Skip to content
Snippets Groups Projects
Commit d9056496 authored by Nico Castro-Folker's avatar Nico Castro-Folker
Browse files

Update Science.hpp to include declaration of compute_lambda2 function

parent 1425e905
No related branches found
No related tags found
No related merge requests found
...@@ -35,6 +35,14 @@ void compute_vorticity(TArrayn::DTArray & vortx, TArrayn::DTArray & vorty, TArra ...@@ -35,6 +35,14 @@ void compute_vorticity(TArrayn::DTArray & vortx, TArrayn::DTArray & vorty, TArra
TArrayn::DTArray & u, TArrayn::DTArray & v, TArrayn::DTArray & w, TArrayn::DTArray & u, TArrayn::DTArray & v, TArrayn::DTArray & w,
TArrayn::Grad * gradient_op, const string * grid_type); TArrayn::Grad * gradient_op, const string * grid_type);
// Baroclinic Vorticity
void compute_baroclinic_vort_x(TArrayn::DTArray & barovortx, TArrayn::DTArray & T,
TArrayn::Grad * gradient_op, const string * grid_type);
void compute_baroclinic_vort_y(TArrayn::DTArray & barovorty, TArrayn::DTArray & T,
TArrayn::Grad * gradient_op, const string * grid_type);
void compute_baroclinic_vort(TArrayn::DTArray & barovort, TArrayn::DTArray & temp,
TArrayn::DTArray & T, TArrayn::Grad * gradient_op, const string * grid_type, bool v_exist);
// Enstrophy density // Enstrophy density
void enstrophy_density(TArrayn::DTArray & enst, TArrayn::DTArray & u, TArrayn::DTArray & v, void enstrophy_density(TArrayn::DTArray & enst, TArrayn::DTArray & u, TArrayn::DTArray & v,
TArrayn::DTArray & w, TArrayn::Grad * gradient_op, const string * grid_type, TArrayn::DTArray & w, TArrayn::Grad * gradient_op, const string * grid_type,
...@@ -122,6 +130,14 @@ void R_invt(TArrayn::DTArray & R, TArrayn::DTArray & u, ...@@ -122,6 +130,14 @@ void R_invt(TArrayn::DTArray & R, TArrayn::DTArray & u,
TArrayn::DTArray & temp2, TArrayn::Grad * gradient_op, TArrayn::DTArray & temp2, TArrayn::Grad * gradient_op,
const string * grid_type, bool v_exist); const string * grid_type, bool v_exist);
// lambda2, second eigenvalue of S^2+Omega^2
void compute_lambda2(TArrayn::DTArray & lambda2, TArrayn::DTArray & u,
TArrayn::DTArray & v, TArrayn::DTArray & w, TArrayn::DTArray & temp1,
TArrayn::DTArray & temp2, TArrayn::Grad * gradient_op,
const string * grid_type, TArrayn::DTArray & A11, TArrayn::DTArray & A12,
TArrayn::DTArray & A13, TArrayn::DTArray & A22, TArrayn::DTArray & A23,
TArrayn::DTArray & A33);
// Equation of state for seawater, polynomial fit from // Equation of state for seawater, polynomial fit from
// Brydon, Sun, Bleck (1999) (JGR) // Brydon, Sun, Bleck (1999) (JGR)
...@@ -146,6 +162,49 @@ inline double eqn_of_state(double T, double S){ ...@@ -146,6 +162,49 @@ inline double eqn_of_state(double T, double S){
// Define a Blitz-friendly operator // Define a Blitz-friendly operator
BZ_DECLARE_FUNCTION2(eqn_of_state) BZ_DECLARE_FUNCTION2(eqn_of_state)
// Derivative of the equation of state WRT T
inline double eqn_of_state_dT(double T){
// Returbs the derivative of the eqn_of_state function above WRT T
// temperature T (degrees celsius), S assumed to be zero
// Constants are from table 4 of the above paper, at pressure 0
// (This is appropriate since this is currently an incompressible
// model)
// Kept the same names for these constants as in eqn_of_state
const double c2 = 5.10768e-2; // T term, constant after derivative
const double c4 = -7.40849e-3; // T^2 term, 2T after derivative
const double c5 = -3.01036e-3; // ST term, S after derivative
const double c6 = 3.32267e-5; // T^3 term, 3t^2 after derivative
const double c7 = 3.21931e-5; // ST^2 term, 2ST after derivative
return c2 + 2*c4*T + 3*c6*T*T;
}
// Define a Blitz-friendly operator
BZ_DECLARE_FUNCTION(eqn_of_state_dT)
// Derivative of the equation of state WRT S
inline double eqn_of_state_dS(double T, double S){
// Returns the density anomaly (kg/m^3) for water at the given
// temperature T (degrees celsius) and salinity S (PSU?)
// Constants are from table 4 of the above paper, at pressure 0
// (This is appropriate since this is currently an incompressible
// model)
const double c3 = 8.05999e-1; // S term, constant after derivative
const double c5 = -3.01036e-3; // ST term, T after derivative
const double c7 = 3.21931e-5; // ST^2 term, T^2 after derivative
return c3 + c5*T + c7*T*T;
}
// Define a Blitz-friendly operator
BZ_DECLARE_FUNCTION2(eqn_of_state_dS)
inline double eqn_of_state_t(double T){ inline double eqn_of_state_t(double T){
// Specialize for freshwater with S=0 // Specialize for freshwater with S=0
return eqn_of_state(T,0.0); return eqn_of_state(T,0.0);
......
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment