


Convert tidal amplitude and phase lag (ap-) parameters into tidal ellipse
(ep-) parameters. Please refer to ep2app for its inverse function.
Usage:
[SEMA, ECC, INC, PHA, w]=ap2ep(Au, PHIu, Av, PHIv, plot_demo)
where:
Au, PHIu, Av, PHIv are the amplitudes and phase lags (in degrees) of
u- and v- tidal current components. They can be vectors or
matrices or multidimensional arrays.
plot_demo is an optional argument, when it is supplied as an array
of indices, say [i j k l], the program will plot an ellipse
corresponding to Au(i,j, k, l), PHIu(i,j,k,l), Av(i,j,k,l), and
PHIv(i,j,k,l);
Any number of dimensions are allowed as long as your computer
resource can handle.
SEMA: Semi-major axes, or the maximum speed;
ECC: Eccentricity, the ratio of semi-minor axis over
the semi-major axis; its negative value indicates that the ellipse
is traversed in clockwise direction.
INC: Inclination, the angles (in degrees) between the semi-major
axes and u-axis.
PHA: Phase angles, the time (in angles and in degrees) when the
tidal currents reach their maximum speeds, (i.e.
PHA=omega*tmax).
These four ep-parameters will have the same dimensionality
(i.e., vectors, or matrices) as the input ap-parameters.
w: Optional. If it is requested, it will be output as matrices
whose rows allow for plotting ellipses and whose columns are
for different ellipses corresponding columnwise to SEMA. For
example, plot(real(w(1,:)), imag(w(1,:))) will let you see
the first ellipse. You may need to use squeeze function when
w is a more than two dimensional array. See example.m.
Document: tidal_ellipse.ps
Revisions: May 2002, by Zhigang Xu, --- adopting Foreman's northern
semi major axis convention.
For a given ellipse, its semi-major axis is undetermined by 180. If we borrow
Foreman's terminology to call a semi major axis whose direction lies in a range of
[0, 180) as the northern semi-major axis and otherwise as a southern semi major
axis, one has freedom to pick up either northern or southern one as the semi major
axis without affecting anything else. Foreman (1977) resolves the ambiguity by
always taking the northern one as the semi-major axis. This revision is made to
adopt Foreman's convention. Note the definition of the phase, PHA, is still
defined as the angle between the initial current vector, but when converted into
the maximum current time, it may not give the time when the maximum current first
happens; it may give the second time that the current reaches the maximum
(obviously, the 1st and 2nd maximum current times are half tidal period apart)
depending on where the initial current vector happen to be and its rotating sense.
Version 2, May 2002

0001 function [SEMA, ECC, INC, PHA, w, TWOCIR]=ap2ep(Au, PHIu, Av, PHIv, plot_demo) 0002 % 0003 % Convert tidal amplitude and phase lag (ap-) parameters into tidal ellipse 0004 % (ep-) parameters. Please refer to ep2app for its inverse function. 0005 % 0006 % Usage: 0007 % 0008 % [SEMA, ECC, INC, PHA, w]=ap2ep(Au, PHIu, Av, PHIv, plot_demo) 0009 % 0010 % where: 0011 % 0012 % Au, PHIu, Av, PHIv are the amplitudes and phase lags (in degrees) of 0013 % u- and v- tidal current components. They can be vectors or 0014 % matrices or multidimensional arrays. 0015 % 0016 % plot_demo is an optional argument, when it is supplied as an array 0017 % of indices, say [i j k l], the program will plot an ellipse 0018 % corresponding to Au(i,j, k, l), PHIu(i,j,k,l), Av(i,j,k,l), and 0019 % PHIv(i,j,k,l); 0020 % 0021 % Any number of dimensions are allowed as long as your computer 0022 % resource can handle. 0023 % 0024 % SEMA: Semi-major axes, or the maximum speed; 0025 % ECC: Eccentricity, the ratio of semi-minor axis over 0026 % the semi-major axis; its negative value indicates that the ellipse 0027 % is traversed in clockwise direction. 0028 % INC: Inclination, the angles (in degrees) between the semi-major 0029 % axes and u-axis. 0030 % PHA: Phase angles, the time (in angles and in degrees) when the 0031 % tidal currents reach their maximum speeds, (i.e. 0032 % PHA=omega*tmax). 0033 % 0034 % These four ep-parameters will have the same dimensionality 0035 % (i.e., vectors, or matrices) as the input ap-parameters. 0036 % 0037 % w: Optional. If it is requested, it will be output as matrices 0038 % whose rows allow for plotting ellipses and whose columns are 0039 % for different ellipses corresponding columnwise to SEMA. For 0040 % example, plot(real(w(1,:)), imag(w(1,:))) will let you see 0041 % the first ellipse. You may need to use squeeze function when 0042 % w is a more than two dimensional array. See example.m. 0043 % 0044 % Document: tidal_ellipse.ps 0045 % 0046 % Revisions: May 2002, by Zhigang Xu, --- adopting Foreman's northern 0047 % semi major axis convention. 0048 % 0049 % For a given ellipse, its semi-major axis is undetermined by 180. If we borrow 0050 % Foreman's terminology to call a semi major axis whose direction lies in a range of 0051 % [0, 180) as the northern semi-major axis and otherwise as a southern semi major 0052 % axis, one has freedom to pick up either northern or southern one as the semi major 0053 % axis without affecting anything else. Foreman (1977) resolves the ambiguity by 0054 % always taking the northern one as the semi-major axis. This revision is made to 0055 % adopt Foreman's convention. Note the definition of the phase, PHA, is still 0056 % defined as the angle between the initial current vector, but when converted into 0057 % the maximum current time, it may not give the time when the maximum current first 0058 % happens; it may give the second time that the current reaches the maximum 0059 % (obviously, the 1st and 2nd maximum current times are half tidal period apart) 0060 % depending on where the initial current vector happen to be and its rotating sense. 0061 % 0062 % Version 2, May 2002 0063 0064 if nargin < 5 0065 plot_demo=0; % by default, no plot for the ellipse 0066 end 0067 0068 0069 % Assume the input phase lags are in degrees and convert them in radians. 0070 PHIu = PHIu/180*pi; 0071 PHIv = PHIv/180*pi; 0072 0073 % Make complex amplitudes for u and v 0074 i = sqrt(-1); 0075 u = Au.*exp(-i*PHIu); 0076 v = Av.*exp(-i*PHIv); 0077 0078 % Calculate complex radius of anticlockwise and clockwise circles: 0079 wp = (u+i*v)/2; % for anticlockwise circles 0080 wm = conj(u-i*v)/2; % for clockwise circles 0081 % and their amplitudes and angles 0082 Wp = abs(wp); 0083 Wm = abs(wm); 0084 THETAp = angle(wp); 0085 THETAm = angle(wm); 0086 0087 % calculate ep-parameters (ellipse parameters) 0088 SEMA = Wp+Wm; % Semi Major Axis, or maximum speed 0089 SEMI = Wp-Wm; % Semin Minor Axis, or minimum speed 0090 ECC = SEMI./SEMA; % Eccentricity 0091 0092 PHA = (THETAm-THETAp)/2; % Phase angle, the time (in angle) when 0093 % the velocity reaches the maximum 0094 INC = (THETAm+THETAp)/2; % Inclination, the angle between the 0095 % semi major axis and x-axis (or u-axis). 0096 0097 % convert to degrees for output 0098 PHA = PHA/pi*180; 0099 INC = INC/pi*180; 0100 THETAp = THETAp/pi*180; 0101 THETAm = THETAm/pi*180; 0102 0103 %map the resultant angles to the range of [0, 360]. 0104 PHA=mod(PHA+360, 360); 0105 INC=mod(INC+360, 360); 0106 0107 % Mar. 2, 2002 Revision by Zhigang Xu (REVISION_1) 0108 % Change the southern major axes to northern major axes to conform the tidal 0109 % analysis convention (cf. Foreman, 1977, p. 13, Manual For Tidal Currents 0110 % Analysis Prediction, available in www.ios.bc.ca/ios/osap/people/foreman.htm) 0111 k = fix(INC/180); 0112 INC = INC-k*180; 0113 PHA = PHA+k*180; 0114 PHA = mod(PHA, 360); 0115 0116 % plot at the request 0117 if nargout >= 5 0118 ndot=36; 0119 dot=2*pi/ndot; 0120 ot=[0:dot:2*pi-dot]; 0121 w=wp(:)*exp(i*ot)+wm(:)*exp(-i*ot); 0122 w=reshape(w, [size(Au) ndot]); 0123 end 0124 0125 if any(plot_demo) 0126 plot_ell(SEMA, ECC, INC, PHA, plot_demo); 0127 end 0128 0129 if nargout == 6 0130 TWOCIR=struct('Wp', Wp, 'THETAp', THETAp, 'wp', ... 0131 wp, 'Wm', Wm, 'THETAm', THETAm, 'wm', wm, 'ot', ot, 'dot', dot); 0132 end 0133 0134 0135 %Authorship Copyright: 0136 % 0137 % The author retains the copyright of this program, while you are welcome 0138 % to use and distribute it as long as you credit the author properly and respect 0139 % the program name itself. Particularly, you are expected to retain the original 0140 % author's name in this original version or any of its modified version that 0141 % you might make. You are also expected not to essentially change the name of 0142 % the programs except for adding possible extension for your own version you 0143 % might create, e.g. ap2ep_xx is acceptable. Any suggestions are welcome and 0144 % enjoy my program(s)! 0145 % 0146 % 0147 %Author Info: 0148 %_______________________________________________________________________ 0149 % Zhigang Xu, Ph.D. 0150 % (pronounced as Tsi Gahng Hsu) 0151 % Research Scientist 0152 % Coastal Circulation 0153 % Bedford Institute of Oceanography 0154 % 1 Challenge Dr. 0155 % P.O. Box 1006 Phone (902) 426-2307 (o) 0156 % Dartmouth, Nova Scotia Fax (902) 426-7827 0157 % CANADA B2Y 4A2 email xuz@dfo-mpo.gc.ca 0158 %_______________________________________________________________________ 0159 % 0160 % Release Date: Nov. 2000, Revised on May. 2002 to adopt Foreman's northern semi 0161 % major axis convention. 0162 0163