XBeach
|
00001 module readwind_module 00002 implicit none 00003 save 00004 private 00005 public readwind 00006 contains 00007 subroutine readwind (s,par) 00008 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00009 ! Copyright (C) 2007 UNESCO-IHE, WL|Delft Hydraulics and Delft University ! 00010 ! Dano Roelvink, Ap van Dongeren, Ad Reniers, Jamie Lescinski, ! 00011 ! Jaap van Thiel de Vries, Robert McCall ! 00012 ! ! 00013 ! d.roelvink@unesco-ihe.org ! 00014 ! UNESCO-IHE Institute for Water Education ! 00015 ! P.O. Box 3015 ! 00016 ! 2601 DA Delft ! 00017 ! The Netherlands ! 00018 ! ! 00019 ! This library is free software; you can redistribute it and/or ! 00020 ! modify it under the terms of the GNU Lesser General Public ! 00021 ! License as published by the Free Software Foundation; either ! 00022 ! version 2.1 of the License, or (at your option) any later version. ! 00023 ! ! 00024 ! This library is distributed in the hope that it will be useful, ! 00025 ! but WITHOUT ANY WARRANTY; without even the implied warranty of ! 00026 ! MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU ! 00027 ! Lesser General Public License for more details. ! 00028 ! ! 00029 ! You should have received a copy of the GNU Lesser General Public ! 00030 ! License along with this library; if not, write to the Free Software ! 00031 ! Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 ! 00032 ! USA ! 00033 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!! 00034 use params, only: parameters 00035 use spaceparams 00036 use xmpi_module 00037 use logging_module, only: writelog, report_file_read_error 00038 use constants, only: pi 00039 00040 implicit none 00041 00042 type(spacepars),target :: s 00043 type(parameters), intent(in) :: par 00044 00045 integer :: i 00046 integer :: nwind,io 00047 real*8,dimension(3) :: temp 00048 00049 00050 00051 if(.not. xmaster) return 00052 00053 io = 0 00054 nwind = 0 00055 00056 allocate(s%windsu(s%nx+1,s%ny+1)) 00057 allocate(s%windnv(s%nx+1,s%ny+1)) 00058 s%windsu=0.d0 00059 s%windnv=0.d0 00060 if (par%windfile==' ') then ! Stationary wind 00061 00062 s%windlen=1 00063 allocate(s%windinpt(1)) 00064 allocate(s%windvelts (1)) 00065 allocate(s%winddirts (1)) 00066 allocate(s%windxts(1)) 00067 allocate(s%windyts(1)) 00068 00069 s%windinpt=par%tstop 00070 s%windvelts=par%windv 00071 s%winddirts=(270.d0-par%windth-s%alfa)*par%px/180.d0 00072 ! Alfa is East to X still ?? 00073 s%windxts = dcos(s%winddirts)*s%windvelts 00074 s%windyts = dsin(s%winddirts)*s%windvelts 00075 ! If stationary wind then we will calculate now what the wind velocity is everywhere now 00076 s%windsu = s%windxts(1)*dcos(s%alfau) + s%windyts(1)*dsin(s%alfau) 00077 s%windnv = s%windyts(1)*dcos(s%alfav-0.5d0*par%px) - s%windxts(1)*dsin(s%alfav-0.5d0*par%px) 00078 else ! Non-stationary wind 00079 call writelog('ls','','readwind: reading wind time series from ',trim(par%windfile),' ...') 00080 open(31,file=par%windfile) 00081 do while (io==0) 00082 nwind=nwind+1 00083 read(31,*,IOSTAT=io) temp 00084 enddo 00085 rewind(31) 00086 s%windlen=nwind-1 00087 allocate(s%windinpt(s%windlen)) 00088 allocate(s%windvelts(s%windlen)) 00089 allocate(s%winddirts(s%windlen)) 00090 allocate(s%windxts(s%windlen)) 00091 allocate(s%windyts(s%windlen)) 00092 do i=1,s%windlen 00093 read(31,*,IOSTAT=io) s%windinpt(i),s%windvelts(i),s%winddirts(i) 00094 if (io .ne. 0) then 00095 call report_file_read_error(par%windfile) 00096 endif 00097 enddo 00098 ! to cartesian radians 00099 s%winddirts=(270.d0-s%winddirts-s%alfa)*par%px/180.d0 00100 s%windxts = dcos(s%winddirts)*s%windvelts 00101 s%windyts = dsin(s%winddirts)*s%windvelts 00102 close(31) 00103 if (par%morfacopt==1) s%windinpt = s%windinpt / max(par%morfac,1.d0) 00104 if (s%windinpt(s%windlen)<par%tstop) then 00105 call writelog('els','','Wind condition time series too short. Stopping calculation') 00106 call halt_program 00107 endif 00108 endif 00109 00110 end subroutine readwind 00111 00112 end module readwind_module