XBeach
C:/repositories/XBeach/trunk/src/xbeachlibrary/readtide.F90
Go to the documentation of this file.
00001 module readtide_module
00002    implicit none
00003    save
00004    private
00005    public readtide
00006 contains
00007    subroutine readtide (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 
00038       use logging_module, only: writelog, report_file_read_error
00039  
00040 
00041       implicit none
00042 
00043       type(spacepars),target                  :: s
00044       type(parameters)                        :: par
00045 
00046       integer                             :: i
00047       integer                             :: io,ntide
00048       real*8,dimension(par%tideloc+1)     :: temp
00049 
00050       ! this must only work on master
00051       if(.not. xmaster) return
00052       if(.true.) then
00053          if (par%tideloc .eq. 0) then
00054             s%tidelen = 2
00055             allocate(s%tideinpt(s%tidelen))
00056             allocate(s%tideinpz(s%tidelen,par%tideloc))
00057             return
00058          endif
00059 
00060          io = 0
00061          ntide = 0
00062 
00063          call writelog('ls','','readtide: reading tide time series from ',trim(par%zs0file),' ...')
00064          open(31,file=par%zs0file)
00065          do while (io==0)
00066             ntide=ntide+1
00067             read(31,*,IOSTAT=io) temp
00068          enddo
00069          rewind(31)
00070 
00071          s%tidelen=ntide-1
00072 
00073          allocate(s%tideinpt(s%tidelen))
00074          allocate(s%tideinpz(s%tidelen,par%tideloc))
00075          s%tideinpz = 0.0d0
00076          do i=1,s%tidelen
00077             read(31,*,iostat=io)s%tideinpt(i),s%tideinpz(i,:)
00078             if (io .ne. 0) then
00079                call report_file_read_error(par%zs0file)
00080             endif
00081          end do
00082          close(31)
00083 
00084          if (par%morfacopt==1) s%tideinpt = s%tideinpt / max(par%morfac,1.d0)
00085          if (s%tideinpt(s%tidelen)<par%tstop) then
00086             call writelog('els','','Tide condition time series too short. Stopping calculation')
00087             call halt_program
00088          endif
00089 
00090       endif
00091 
00092    end subroutine readtide
00093 
00094 end module readtide_module
 All Classes Files Functions Variables Defines