Home > wafo > wavemodels > jhvpdf2.m

jhvpdf2

PURPOSE ^

Joint (Vcf,Hd) PDF for linear waves with a JONSWAP spectrum.

SYNOPSIS ^

[f,Hrms,Vrms,fA,fB] = jhvpdf2(Hd,Vcf,Hm0,Tp,gam,normalizedInput)

DESCRIPTION ^

 JHVPDF2 Joint (Vcf,Hd) PDF for linear waves with a JONSWAP spectrum. 
  
   CALL: f = jhvpdf2(Hd,Vcf,Hm0,Tp,gamma) 
   
   f     = pdf struct evaluated at meshgrid(Vcf,Hd) 
   Hd    = zero down crossing wave height 
   Vcf   = crest front velocity 
   Hm0   = significant wave height 
   Tp    = Spectral peak period  
   Gamma = Peakedness parameter of the JONSWAP spectrum 
  
  JHVPDF2 approximates the joint distribution of (Vcf, Hd), i.e., crest 
  front velocity (Ac/Tcf) and wave height, for a Gaussian process with a 
  JONSWAP spectral density. The empirical parameters of the model is 
  fitted by least squares to simulated (Vcf,Hd) data for 13 classes of 
  GAMMA between 1 and 7. About 100000 zero-downcrossing waves 
  were simulated for each class. 
  JHVPDF2 is restricted to the following range for GAMMA:  
   1 <= GAMMA <= 7  
  
  Example: 
  Hm0 = 6;Tp = 9; gam=3.5 
  h = linspace(0,4*Hm0/sqrt(2))';  
  v = linspace(0,4*2*Hm0/Tp)'; 
  f = jhvpdf2(h,v,Hm0,Tp,gam); 
  w = linspace(0,40,5*1024+1).'; 
  S = jonswap(w,[Hm0, Tp, gam]); 
  dt = .3; 
  x = spec2sdat(S,80000,dt); rate = 4; 
  [vi,hi] = dat2steep(x,rate,1); 
  fk = kdebin([vi,hi],'epan',[],[],.5,128);  
   fk.title = f.title; fk.labx = f.labx;   
  plot(vi,hi,'.'), hold on 
  pdfplot(f),pdfplot(fk,'r'),hold off 
  
  See also  thvpdf

CROSS-REFERENCE INFORMATION ^

This function calls: This function is called by:

SOURCE CODE ^

001 function [f,Hrms,Vrms,fA,fB] = jhvpdf2(Hd,Vcf,Hm0,Tp,gam,normalizedInput) 
002 %JHVPDF2 Joint (Vcf,Hd) PDF for linear waves with a JONSWAP spectrum. 
003 % 
004 %  CALL: f = jhvpdf2(Hd,Vcf,Hm0,Tp,gamma) 
005 %  
006 %  f     = pdf struct evaluated at meshgrid(Vcf,Hd) 
007 %  Hd    = zero down crossing wave height 
008 %  Vcf   = crest front velocity 
009 %  Hm0   = significant wave height 
010 %  Tp    = Spectral peak period  
011 %  Gamma = Peakedness parameter of the JONSWAP spectrum 
012 % 
013 % JHVPDF2 approximates the joint distribution of (Vcf, Hd), i.e., crest 
014 % front velocity (Ac/Tcf) and wave height, for a Gaussian process with a 
015 % JONSWAP spectral density. The empirical parameters of the model is 
016 % fitted by least squares to simulated (Vcf,Hd) data for 13 classes of 
017 % GAMMA between 1 and 7. About 100000 zero-downcrossing waves 
018 % were simulated for each class. 
019 % JHVPDF2 is restricted to the following range for GAMMA:  
020 %  1 <= GAMMA <= 7  
021 % 
022 % Example: 
023 % Hm0 = 6;Tp = 9; gam=3.5 
024 % h = linspace(0,4*Hm0/sqrt(2))';  
025 % v = linspace(0,4*2*Hm0/Tp)'; 
026 % f = jhvpdf2(h,v,Hm0,Tp,gam); 
027 % w = linspace(0,40,5*1024+1).'; 
028 % S = jonswap(w,[Hm0, Tp, gam]); 
029 % dt = .3; 
030 % x = spec2sdat(S,80000,dt); rate = 4; 
031 % [vi,hi] = dat2steep(x,rate,1); 
032 % fk = kdebin([vi,hi],'epan',[],[],.5,128);  
033 %  fk.title = f.title; fk.labx = f.labx;   
034 % plot(vi,hi,'.'), hold on 
035 % pdfplot(f),pdfplot(fk,'r'),hold off 
036 % 
037 % See also  thvpdf 
038  
039 % Reference   
040 % P. A. Brodtkorb (2004),   
041 % The Probability of Occurrence of Dangerous Wave Situations at Sea. 
042 % Dr.Ing thesis, Norwegian University of Science and Technolgy, NTNU, 
043 % Trondheim, Norway.    
044    
045 error(nargchk(3,6,nargin)) 
046  
047 if (nargin < 6|isempty(normalizedInput)),  normalizedInput  = 0;end 
048 if (nargin < 4|isempty(Tp)),  Tp  = 8;end 
049 if (nargin < 3|isempty(Hm0)), Hm0 = 6;end 
050 if (nargin < 5|isempty(gam)) 
051    gam = getjonswappeakedness(Hm0,Tp); 
052 end 
053 displayWarning = 1; 
054 if displayWarning 
055   if any(any(Tp>5*sqrt(Hm0) | Tp<3.6*sqrt(Hm0))) 
056     disp('Warning: Hm0,Tp is outside the JONSWAP range') 
057     disp('The validity of the parameters returned are questionable') 
058   end 
059 end 
060  
061 [V,H] = meshgrid(Vcf,Hd); 
062  
063 f = createpdf(2); 
064 [f.f,Hrms,Vrms,varargout{1:nargout-1}]  = jhvpdf(H,V,Hm0,Tp,gam,normalizedInput); 
065 f.x = {Vcf(:),Hd(:)}; 
066   
067 if (normalizedInput) 
068   f.labx = {'Vcf', 'Hd'}; 
069   f.norm = 1; 
070 else 
071   f.norm = 0; 
072   f.labx = {'Vcf [m/s]', 'Hd [m]'}; 
073 end 
074 f.title = 'Joint distribution of (Hd,Vcf) in time'; 
075 f.note = ['Jonswap Hm0=' num2str(Hm0) ' Tp = ' num2str(Tp) ' Gamma = ' num2str(gam)]; 
076 [f.cl,f.pl] = qlevels(f.f); 
077 return 
078

Mathematical Statistics
Centre for Mathematical Sciences
Lund University with Lund Institute of Technology

Comments or corrections to the WAFO group


Generated on Thu 06-Oct-2005 02:21:16 for WAFO by m2html © 2003