#include <boost/test/unit_test.hpp>
#include <usml/waveq3d/waveq3d.h>
#include <usml/types/seq_vector.h>
#include <iostream>
#include <iomanip>
#include <fstream>
BOOST_AUTO_TEST_SUITE(refraction_test)
using namespace boost::unit_test;
using namespace usml::waveq3d;
static const double time_step = 0.1 ;
static const seq_log freq(10e3, 10e3, 1);
BOOST_AUTO_TEST_CASE(refraction_isovelocity) {
cout << "=== refraction_test: refraction_isovelocity ===" << endl;
const double c0 = 1500.0;
profile_model* profile = new profile_linear(c0);
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(5000.0);
ocean_model ocean(surface, bottom, profile);
wposition1 pos(45.0, -45.0, -1000.0);
seq_linear de(0.0, 1.0, 1);
seq_linear az(0.0, 30.0, 90.0);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << " Hz" << endl;
double rho = pos.rho();
double st = sin(pos.theta());
double ct = cos(pos.theta());
double sp = sin(pos.phi());
double cp = cos(pos.phi());
const double x0 = rho * st * cp;
const double y0 = rho * st * sp;
const double z0 = rho * ct;
vector<double> east = az;
east = sin(to_radians(east));
vector<double> north = az;
north = cos(to_radians(north));
vector<double> x_dir = -sp * east - ct * cp * north;
vector<double> y_dir = cp * east - ct * sp * north;
vector<double> z_dir = st * north;
const double r0 = rho;
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_isovelocity.csv";
std::ofstream os(name);
cout << "writing tables to " << name << endl;
os << "time,x,y,z,d";
for (size_t n = 0; n < az.size(); ++n) {
os << ",x" << az(n)
<< ",y" << az(n)
<< ",z" << az(n)
<< ",d" << az(n);
}
os << endl;
double max_error = 0.0;
while (wave.curr()->position.altitude(0,0) < -10.0) {
double time = wave.time();
vector<double> x;
x = x0 + c0 * time * x_dir;
vector<double> y;
y = y0 + c0 * time * y_dir;
vector<double> z;
z = z0 + c0 * time * z_dir;
double d_analytic = r0 / cos(atan(c0 * time / r0))
- wposition::earth_radius;
os << time << "," << x(0) << "," << y(0) << "," << z(0)
<< "," << d_analytic;
for (size_t n = 0; n < az.size(); ++n) {
rho = wave.curr()->position.rho(0, n);
st = sin(wave.curr()->position.theta(0, n));
ct = cos(wave.curr()->position.theta(0, n));
sp = sin(wave.curr()->position.phi(0, n));
cp = cos(wave.curr()->position.phi(0, n));
double x_model = rho * st * cp;
double y_model = rho * st * sp;
double z_model = rho * ct;
double dx = x_model - x(n);
max_error = max(max_error, fabs(dx));
double dy = y_model - y(n);
max_error = max(max_error, fabs(dy));
double dz = z_model - z(n);
max_error = max(max_error, fabs(dz));
os << "," << x_model << "," << y_model << "," << z_model;
BOOST_CHECK(fabs(dx) < 1e-3);
BOOST_CHECK(fabs(dy) < 1e-3);
BOOST_CHECK(fabs(dz) < 1e-3);
double alt_model = wave.curr()->position.altitude(0, n);
double dd = alt_model - d_analytic;
max_error = max(max_error, fabs(dd));
os << "," << alt_model;
BOOST_CHECK(fabs(dd) < 1e-3);
}
os << endl;
wave.step();
}
cout << "wave breaks surface around " << wave.time() << " secs" << endl
<< "max error = " << max_error << " meters" << endl;
}
BOOST_AUTO_TEST_CASE(refraction_great_circle) {
cout << "=== refraction_test: refraction_great_circle ===" << endl;
const double c0 = 1500.0;
profile_model* profile = new profile_linear(c0);
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(5000.0);
ocean_model ocean(surface, bottom, profile);
profile->flat_earth(true);
double lat1 = 45.0;
double lng1 = -45.0;
double alt1 = -1000.0;
wposition1 pos(lat1, lng1, alt1);
seq_linear de(0.0, 1.0, 1);
seq_linear az(0.0, 30.0, 90.0);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << " Hz" << endl;
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_great_circle.csv";
std::ofstream os(name);
cout << "writing tables to " << name << endl;
os << "time";
for (size_t n = 0; n < az.size(); ++n) {
os << ",lat" << az(n) << ",long" << az(n)
<< ",dalt" << az(n) << ",dbear" << az(n);
}
os << endl;
lat1 *= M_PI / 180;
lng1 *= M_PI / 180;
double max_d_alt = 0.0;
double max_d_tc1 = 0.0;
while (wave.time() < 1000.0) {
wave.step();
os << wave.time();
for (size_t n = 0; n < az.size(); ++n) {
double alt2 = wave.curr()->position.altitude(0, n);
double d_alt = alt2 - alt1;
max_d_alt = max(max_d_alt, fabs(d_alt));
double lat2 = wave.curr()->position.latitude(0, n)
* M_PI / 180;
double lng2 = wave.curr()->position.longitude(0, n)
* M_PI / 180;
double tc1 = -atan2(sin(lng1 - lng2) * cos(lat2),
cos(lat1) * sin(lat2) - sin(lat1) * cos(lat2) * cos(lng1 - lng2))
* 180 / M_PI;
double d_tc1 = tc1 - az(n);
max_d_tc1 = max(max_d_tc1, fabs(d_tc1));
os << "," << wave.curr()->position.latitude(0, n)
<< "," << wave.curr()->position.longitude(0, n)
<< "," << d_alt << "," << d_tc1;
BOOST_CHECK(fabs(d_alt) < 1e-3);
BOOST_CHECK(fabs(d_tc1) < 1e-3);
}
os << endl;
}
cout << "wave propagates for " << wave.time() << " secs" << endl
<< "max error = " << max_d_alt << " meters and "
<< max_d_tc1 << " degrees" << endl;
}
BOOST_AUTO_TEST_CASE(refraction_linear) {
cout << "=== refraction_test: refraction_linear ===" << endl;
const double angle = 0.0;
const double ang = to_radians( angle );
const double z0 = 1000.0;
const double c0 = 1500.0;
const double g0 = 0.016;
const double a0 = cos(ang) / (c0 + g0 * z0);
const double sinT = sin(ang);
profile_model* profile = new profile_linear(c0, g0);
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(5000.0);
ocean_model ocean(surface, bottom, profile);
profile->flat_earth(true);
wposition1 pos(0.0, 0.0, -z0);
seq_linear de(angle, 0.0, 1);
seq_linear az(90.0, 0.0, 1);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << " Hz" << endl;
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_linear.csv";
std::ofstream os(name);
cout << "writing tables to " << name << endl;
os << "t," << "r," << "theta," << "phi," << "z," << "rng(m),"
<< "rd," << "thd," << "phid,"
<< "alpha," << "beta," << "gamma,"
<< "alphad," << "betad," << "gammad,"
<< "c," << "dcdz,"
<< endl;
os << std::scientific << std::showpoint << std::setprecision(18);
wposition prev(1, 1);
prev.rho(0, 0, pos.rho());
prev.theta(0, 0, pos.theta());
prev.phi(0, 0, pos.phi());
double max_error = 0.0;
double r = 0.0;
double z = -z0;
while (wave.time() < 9.0) {
r += 0.5 * (wave.curr()->position.rho(0,0) + prev.rho(0,0))
* (wave.curr()->position.phi(0,0) - prev.phi(0,0));
double agr = a0 * g0 * r + sinT;
z = -(sqrt(1.0 - agr * agr) / a0 - c0) / g0;
prev = wave.curr()->position;
os << wave.time() << ','
<< wave.curr()->position.rho(0,0) << ','
<< wave.curr()->position.theta(0,0) << ','
<< wave.curr()->position.phi(0,0) << ','
<< z << ',' << r << ','
<< wave.curr()->pos_gradient.rho(0,0) << ','
<< wave.curr()->pos_gradient.theta(0,0) << ','
<< wave.curr()->pos_gradient.phi(0,0) << ','
<< wave.curr()->ndirection.rho(0,0) << ','
<< wave.curr()->ndirection.theta(0,0) << ','
<< wave.curr()->ndirection.phi(0,0) << ','
<< wave.curr()->ndir_gradient.rho(0,0) << ','
<< wave.curr()->ndir_gradient.theta(0,0) << ','
<< wave.curr()->ndir_gradient.phi(0,0) << ','
<< wave.curr()->sound_speed(0,0) << ','
<< wave.curr()->sound_gradient.rho(0,0) << endl;
double d_alt = fabs(wave.curr()->position.altitude(0,0) - z);
max_error = max(max_error, d_alt);
BOOST_CHECK(d_alt < 1e-3);
wave.step();
}
cout << "max error = " << max_error << " meters" << endl;
}
BOOST_AUTO_TEST_CASE(refraction_n2_linear) {
cout << "=== refraction_test: refraction_n2_linear ===" << endl;
const double angle = 50.0;
const double ang = to_radians( angle ) ;
const double z0 = 1000.0;
const double c0 = 1550.0;
const double g0 = 1.2 ;
const double a0 = cos(ang) / (c0 / sqrt(1.0 + 2*g0/c0 * z0));
profile_model* profile = new profile_n2(c0, g0) ;
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(5000.0);
ocean_model ocean(surface, bottom, profile);
profile->flat_earth(true);
wposition1 pos(0.0, 0.0, -z0);
seq_linear de(angle, 0.0, 1);
seq_linear az(90.0, 0.0, 1);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << endl;
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_n2_linear.csv";
std::ofstream os(name);
cout << "writing tables to " << name << endl;
os << "t," << "r," << "theta," << "phi," << "z," << "rng(m),"
<< "rd," << "thd," << "phid,"
<< "alpha," << "beta," << "gamma,"
<< "alphad," << "betad," << "gammad,"
<< "c," << "dcdz,"
<< endl;
os << std::scientific << std::showpoint << std::setprecision(18);
wposition prev(1, 1);
prev.rho(0, 0, pos.rho());
prev.theta(0, 0, pos.theta());
prev.phi(0, 0, pos.phi());
double max_error = 0.0;
double r = 0.0;
double z = 0.0;
while (z > -2000.0) {
double time = wave.time();
r += 0.5 * (wave.curr()->position.rho(0,0) + prev.rho(0,0))
* (wave.curr()->position.phi(0,0) - prev.phi(0,0));
z = -(2*g0/c0 / (4.0 * a0 * a0 * c0 * c0) * r * r
- r * tan(ang) + z0);
prev = wave.curr()->position;
os << time << ','
<< wave.curr()->position.rho(0,0) << ','
<< wave.curr()->position.theta(0,0) << ','
<< wave.curr()->position.phi(0,0) << ','
<< z << ',' << r << ','
<< wave.curr()->pos_gradient.rho(0,0) << ','
<< wave.curr()->pos_gradient.theta(0,0) << ','
<< wave.curr()->pos_gradient.phi(0,0) << ','
<< wave.curr()->ndirection.rho(0,0) << ','
<< wave.curr()->ndirection.theta(0,0) << ','
<< wave.curr()->ndirection.phi(0,0) << ','
<< wave.curr()->ndir_gradient.rho(0,0) << ','
<< wave.curr()->ndir_gradient.theta(0,0) << ','
<< wave.curr()->ndir_gradient.phi(0,0) << ','
<< wave.curr()->sound_speed(0,0) << ','
<< wave.curr()->sound_gradient.rho(0,0) << endl;
double d_alt = fabs(wave.curr()->position.altitude(0,0) - z);
max_error = max(max_error, d_alt);
BOOST_CHECK_SMALL( wave.curr()->position.altitude(0,0)-z, 2.0 ) ;
wave.step();
}
cout << "wave reaches " << z << " m depth around "
<< wave.time() << " secs" << endl
<< "max error = " << max_error << " meters" << endl;
}
BOOST_AUTO_TEST_CASE(refraction_catenary) {
cout << "=== refraction_test: refraction_catenary ===" << endl;
const double angle = -3.0;
const double ang = to_radians( angle );
const double z0 = 1000.0;
const double c1 = 1500.0;
const double g1 = 1500.0;
const double sinT2 = sin(ang) * sin(ang);
const double cosT2 = cos(ang) * cos(ang);
profile_model* profile = new profile_catenary(c1, g1, z0);
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(5000.0);
ocean_model ocean(surface, bottom, profile);
profile->flat_earth(true);
wposition1 pos(0.0, 0.0, -z0);
seq_linear de(angle, 0.0, 1);
seq_linear az(90.0, 0.0, 1);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << " Hz" << endl;
wposition prev(1, 1);
prev.rho(0, 0, pos.rho());
prev.theta(0, 0, pos.theta());
prev.phi(0, 0, pos.phi());
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_catenary.csv";
std::ofstream os(name);
cout << "writing tables to " << name << endl;
os << "t," << "r," << "theta," << "phi," << "z," << "rng(m),"
<< "rd," << "thd," << "phid,"
<< "alpha," << "beta," << "gamma,"
<< "alphad," << "betad," << "gammad,"
<< "c," << "dcdz,"
<< endl;
os << std::scientific << std::showpoint << std::setprecision(18);
double max_error = 0.0;
double r = 0.0;
int n = 1;
double sign = -1.0;
while (r < 100000) {
double time = wave.time();
r += 0.5 * (wave.curr()->position.rho(0,0) + prev.rho(0,0))
* (wave.curr()->position.phi(0,0) - prev.phi(0,0));
double t = 2.0 * r / c1;
if (t >= TWO_PI * n) {
++n;
sign *= -1.0;
}
double s = -0.5 * ((sinT2 * cos(t) - 1.0) / cosT2 - 1.0);
double z = -(z0 - sign * c1 * acosh(sqrt(s)));
prev = wave.curr()->position;
os << time << ','
<< wave.curr()->position.rho(0,0) << ','
<< wave.curr()->position.theta(0,0) << ','
<< wave.curr()->position.phi(0,0) << ','
<< z << ',' << r << ','
<< wave.curr()->pos_gradient.rho(0,0) << ','
<< wave.curr()->pos_gradient.theta(0,0) << ','
<< wave.curr()->pos_gradient.phi(0,0) << ','
<< wave.curr()->ndirection.rho(0,0) << ','
<< wave.curr()->ndirection.theta(0,0) << ','
<< wave.curr()->ndirection.phi(0,0) << ','
<< wave.curr()->ndir_gradient.rho(0,0) << ','
<< wave.curr()->ndir_gradient.theta(0,0) << ','
<< wave.curr()->ndir_gradient.phi(0,0) << ','
<< wave.curr()->sound_speed(0,0) << ','
<< wave.curr()->sound_gradient.rho(0,0) << endl;
double d_alt = fabs(wave.curr()->position.altitude(0,0) - z);
max_error = max(max_error, d_alt);
BOOST_CHECK_SMALL( d_alt, 2.0 );
wave.step();
}
cout << "wave propagates for " << wave.time() << " secs" << endl
<< "max error = " << max_error << " meters" << endl;
}
BOOST_AUTO_TEST_CASE(refraction_munk_range) {
cout << "=== refraction_test: refraction_munk_range ===" << endl;
const char* ncname_wave = USML_TEST_DIR "/waveq3d/test/refraction_munk_range.nc";
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_munk_range.csv";
static const double cycle_ranges[] = {
64977.771509,
62686.699943,
60536.790347,
58539.834823,
56706.277890,
55044.418981,
53559.948084,
52255.876772,
51132.827760,
50189.572079,
49423.683193,
48832.195747,
48412.185973,
48161.238557,
48077.771909,
48161.238557,
48412.185973,
48832.195747,
49423.683193,
50189.572079,
51132.827760,
52255.876772,
53559.948084,
55044.418981,
56706.277890,
58539.834823,
60536.790347,
62686.699943,
64977.771509
} ;
profile_munk* profile = new profile_munk();
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(1e4);
ocean_model ocean(surface, bottom, profile);
profile->flat_earth(true);
double lat1 = 45.0;
double lng1 = -45.0;
double alt1 = -1000.0 ;
wposition1 pos(lat1, lng1, alt1);
seq_linear de(-14.0,1.0,14.0);
seq_linear az(0.0, 0.0, 1);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << " Hz" << endl;
std::ofstream os(name);
cout << "writing error data to " << name << endl;
os << "t," << "de," << "Rtheory," << "Rtheory," << "diff," << "diff%,"
<< endl;
os << std::scientific << std::showpoint << std::setprecision(18);
cout << "writing wavefronts to " << ncname_wave << endl;
wave.init_netcdf(ncname_wave);
wave.save_netcdf();
matrix<std::size_t> up( de.size(), 1 ) ;
matrix<std::size_t> low( de.size(), 1 ) ;
up.clear() ;
low.clear() ;
vector<double> loop( de.size() ) ;
loop.clear() ;
double max_error = 0.0 ;
while (wave.time() < 95.0) {
wave.step();
wave.save_netcdf();
for ( size_t d=0 ; d < de.size() ; ++d ) {
const double Hprev = wave.prev()->position.altitude(d,0) - alt1 ;
const double Hcurr = wave.curr()->position.altitude(d,0) - alt1 ;
const double Hnext = wave.next()->position.altitude(d,0) - alt1 ;
if ( Hcurr*Hnext < 0.0 &&
wave.curr()->ndirection.rho(d,0) * de(d) > 0.0 )
{
loop(d) += 1.0 ;
const double Rtheory = loop(d) * cycle_ranges[d] ;
const double Rprev = wave.prev()->position.rho(d,0)
* to_radians(wave.prev()->position.latitude(d,0)-lat1) ;
const double Rcurr = wave.curr()->position.rho(d,0)
* to_radians(wave.curr()->position.latitude(d,0)-lat1) ;
const double Rnext = wave.next()->position.rho(d,0)
* to_radians(wave.next()->position.latitude(d,0)-lat1) ;
double h = (Hnext-Hprev) ;
const double slope = (Rnext-Rprev) / h ;
h *= 0.5 ;
const double curve = (Rnext-2.0*Rcurr+Rprev) / (h*h) ;
const double dx = -Hcurr ;
const double Rmodel = Rcurr + slope*dx + 0.5*curve*dx*dx ;
os << wave.time()
<< "," << de(d) << "," << Rtheory << "," << Rmodel
<< "," << (Rmodel-Rtheory)
<< "," << ((Rmodel-Rtheory)/Rtheory*100.0)
<< endl ;
BOOST_CHECK_CLOSE( Rtheory, Rmodel, 0.01 ) ;
std::size_t dup = wave.curr()->upper(d,0) - up(d,0) ;
std::size_t dlow = wave.curr()->lower(d,0) - low(d,0) ;
BOOST_CHECK_EQUAL( dup, 1 ) ;
BOOST_CHECK_EQUAL( dlow, 1 ) ;
up(d,0) = wave.curr()->upper(d,0) ;
low(d,0) = wave.curr()->lower(d,0) ;
max_error = max( max_error, abs(Rmodel-Rtheory) ) ;
}
}
}
wave.close_netcdf();
cout << "max error = " << max_error << " m" << endl ;
}
BOOST_AUTO_TEST_CASE(refraction_pedersen_range) {
cout << "=== refraction_test: refraction_pedersen_range ===" << endl;
const char* ncname_wave = USML_TEST_DIR "/waveq3d/test/refraction_pedersen_range.nc";
const char* name = USML_TEST_DIR "/waveq3d/test/refraction_pedersen_range.csv";
static const double cycle_ranges[] = {
2115.799965,
2286.528610,
2446.115431,
2593.782977,
2728.811865,
2850.544337,
2958.387468,
3051.815999,
3130.374953,
3193.681828,
3241.428421,
3273.382397,
3289.388370,
3289.368623,
3273.323545,
3241.331594,
} ;
const double c0 = 1550 ;
const double g0 = 1.2 ;
profile_model* profile = new profile_n2( c0, g0 );
boundary_model* surface = new boundary_flat();
boundary_model* bottom = new boundary_flat(1e4);
ocean_model ocean(surface, bottom, profile);
profile->flat_earth(true);
double lat1 = 45.0;
double lng1 = -45.0;
double alt1 = -1000.0 ;
wposition1 pos(lat1, lng1, alt1);
seq_linear de(20.0,2.0,50.0);
seq_linear az(0.0, 0.0, 1);
wave_queue wave(ocean, freq, pos, de, az, time_step);
cout << "time step = " << time_step << " secs"
<< " freq = " << freq(0) << " Hz" << endl;
std::ofstream os(name);
cout << "writing error data to " << name << endl;
os << "t," << "de," << "Rtheory," << "Rtheory," << "diff," << "diff%,"
<< endl;
os << std::scientific << std::showpoint << std::setprecision(18);
cout << "writing wavefronts to " << ncname_wave << endl;
wave.init_netcdf(ncname_wave);
wave.save_netcdf();
double max_error = 0.0 ;
while (wave.time() < 4.0) {
wave.step();
wave.save_netcdf();
for ( size_t d=0 ; d < de.size() ; ++d ) {
const double Hprev = wave.prev()->position.altitude(d,0) - alt1 ;
const double Hcurr = wave.curr()->position.altitude(d,0) - alt1 ;
const double Hnext = wave.next()->position.altitude(d,0) - alt1 ;
if ( Hcurr*Hnext < 0.0 ) {
const double Rtheory = cycle_ranges[d] ;
const double Rprev = wave.prev()->position.rho(d,0)
* to_radians(wave.prev()->position.latitude(d,0)-lat1) ;
const double Rcurr = wave.curr()->position.rho(d,0)
* to_radians(wave.curr()->position.latitude(d,0)-lat1) ;
const double Rnext = wave.next()->position.rho(d,0)
* to_radians(wave.next()->position.latitude(d,0)-lat1) ;
double h = (Hnext-Hprev) ;
const double slope = (Rnext-Rprev) / h ;
h *= 0.5 ;
const double curve = (Rnext-2.0*Rcurr+Rprev) / (h*h) ;
const double dx = -Hcurr ;
const double Rmodel = Rcurr + slope*dx + 0.5*curve*dx*dx ;
os << wave.time()
<< "," << de(d) << "," << Rtheory << "," << Rmodel
<< "," << (Rmodel-Rtheory)
<< "," << ((Rmodel-Rtheory)/Rtheory*100.0)
<< endl ;
BOOST_CHECK_CLOSE( Rtheory, Rmodel, 0.12 );
max_error = max( max_error, abs(Rmodel-Rtheory) ) ;
}
}
}
wave.close_netcdf();
cout << "max error = " << max_error << " m" << endl ;
}
BOOST_AUTO_TEST_CASE( surface_duct_test ) {
cout << "=== refraction_test: surface_duct_test ===" << endl;
const char* ncname_wave = USML_TEST_DIR "/waveq3d/test/refraction_surface_duct.nc";
const char* csvname = USML_TEST_DIR "/waveq3d/test/refraction_surface_duct.csv";
data_grid<double,1>* sound_profile;
seq_vector* axis[1];
axis[0] = new seq_linear( wposition::earth_radius, -0.5, 1000 ) ;
sound_profile = new data_grid<double,1>(axis) ;
size_t index[1];
for(int i=0; i < axis[0]->size(); ++i){
index[0] = i;
double value = 1500.0;
if( (*axis[0])[i] > wposition::earth_radius-150.0 ) { value+=(0.016*i) ;}
if( (*axis[0])[i] <= wposition::earth_radius-150.0 && (*axis[0])[i] > wposition::earth_radius-250.0 ) { value-=(0.1*(i-300)-4.8) ;}
if( (*axis[0])[i] <= wposition::earth_radius-250.0 ) { value-=(0.01*(i-500) + 15.2) ; }
sound_profile->data(index, value);
}
sound_profile->interp_type(0, GRID_INTERP_LINEAR);
sound_profile->edge_limit(0, true);
profile_model* profile = new profile_grid<double,1>(sound_profile);
cout << "writing sound speed profile to " << csvname << endl;
matrix<double> speed(1,1);
wposition test( 1, 1, 45.0, -45.0, 0.0 );
std::ofstream file(csvname);
file << "depth,speed,interp" << endl;
for(int j=0; j <axis[0]->size(); ++j){
index[0] = j;
test.rho( 0, 0, (*axis[0])(j) );
profile->sound_speed( test, &speed );
file << (*axis[0])(j)-wposition::earth_radius << "," << sound_profile->data(index) << "," << speed(0,0) << endl;
}
profile->flat_earth(true);
boundary_model* bottom = new boundary_flat(1e4);
boundary_model* surface = new boundary_flat();
ocean_model ocean(surface, bottom, profile);
double lat = 45.0;
double lon = -45.0;
wposition1 source( lat, lon, -40.0 );
seq_rayfan de ;
seq_linear az( 0.0, 0.0, 1 );
wave_queue wave(ocean, freq, source, de, az, time_step);
cout << "writing wavefronts to " << ncname_wave << endl;
wave.init_netcdf(ncname_wave);
wave.save_netcdf();
while (wave.time() < 25.0) {
wave.step();
wave.save_netcdf();
}
wave.close_netcdf();
delete axis[0] ;
}
BOOST_AUTO_TEST_SUITE_END()