Basilisk CFD
Adaptive Cartesian mesh PDE framework
Loading...
Searching...
No Matches
equilibrium_tide.h File Reference
#include <time.h>
Include dependency graph for equilibrium_tide.h:

Go to the source code of this file.

Data Structures

struct  Etide_t
 

Functions

void equilibrium_tide_constituents (const char *date="2000-01-01 12:00:00", bool Lnodal=false)
 
void equilibrium_tide (scalar tide, double t)
 

Variables

struct { 
 
   Etide_t   Q1 
 
   Etide_t   O1 
 
   Etide_t   K1 
 
   Etide_t   N2 
 
   Etide_t   M2 
 
   Etide_t   S2 
 
   Etide_t   K2 
 
Etide = {0} 
 
const double deg2rad = pi/180.
 

Function Documentation

◆ equilibrium_tide()

void equilibrium_tide ( scalar  tide,
double  t 
)

Computes surface elevation associated with the equilibrium tide

Assumes that the coordinate system is longitude in degrees (x) and latitude in degrees (y).

Note that the tidal constituents must first be defined by calling equilibrium_tide_constituents().

The time t given (in seconds) is relative to the date/time specified when calling equilibrium_tide_constituents().

Compute astronomical equilibrium tide (m) with diurnal and semidiurnal constituents. The long-period constituents and high-order harmonics are neglected.

Definition at line 266 of file equilibrium_tide.h.

References _i, deg2rad, Etide, LINENO, sq(), t, x, and y.

Here is the call graph for this function:

◆ equilibrium_tide_constituents()

void equilibrium_tide_constituents ( const char date = "2000-01-01 12:00:00",
bool  Lnodal = false 
)

Computation of the tidal constituents

The tidal constituents are computed relative to the optional starting date given as argument.

If Lnodal is true then lunar nodal correction is applied.

Compute fundamental astronomical parameters. Adapted from Egbert and Ray, 2017, Table 1.

Set Egbert and Ray time reference for their astronomical parameters (Table 1): days since 2000-01-01:12:00:00

Terrestial time (in centuries) since tide reference date number. Recall that the length of a year is 365.2425 days for the now called Gregorian Calendar (corrected after 15 October 1582).

We also assume that a day is 86400 seconds (need to check this). Note that the original code in ROMS seems to round T to the nearest day??, which we do not do here.*

Compute the harmonic constituents of the equilibrium tide at Greenwich. Adapted from Doodson and Warburg, 1941, Table 1.

Nodal factors "f" and "nu" account for the slow modulation of the tidal constituents due (principally) to the 18.6-year lunar nodal cycle.

Mean longitude of lunar node (Period = 18.6 years).

Mean longitude of the moon (Period = tropical month).

Mean longitude of the sun (Period = tropical year).

Mean longitude of lunar perigee (Period = 8.85 years).

Compute tidal constituent phase "chi" (degrees) at Greenwich meridian.

Set tidal constituent frequency (1/s).

Set tidal constituent amplitude (m).

Set tidal constituent Love number factors (1 + k2 - h2). The Love number is defined as the ratio of the body tide to the height of the static equilibrium tide.

Compute product of amp*f*love.

Definition at line 100 of file equilibrium_tide.h.

References assert, deg2rad, Etide, h, LINENO, N, p, s, T, and x.

Variable Documentation

◆ deg2rad

const double deg2rad = pi/180.

Definition at line 88 of file equilibrium_tide.h.

Referenced by equilibrium_tide(), and equilibrium_tide_constituents().

◆ [struct]

struct { ... } Etide

◆ K1

Etide_t K1

Definition at line 81 of file equilibrium_tide.h.

◆ K2

Etide_t K2

Definition at line 85 of file equilibrium_tide.h.

◆ M2

Etide_t M2

Definition at line 83 of file equilibrium_tide.h.

◆ N2

Etide_t N2

Definition at line 82 of file equilibrium_tide.h.

◆ O1

Etide_t O1

Definition at line 80 of file equilibrium_tide.h.

◆ Q1

Etide_t Q1

Definition at line 79 of file equilibrium_tide.h.

◆ S2

Etide_t S2

Definition at line 84 of file equilibrium_tide.h.