RSTT  3.1.0
Regional Seismic Travel Time
slbm_F_shell.h File Reference
#include "SLBMGlobals.h"

Go to the source code of this file.

Macros

#define SLBM_FORT_SHELL_EXP_IMP
 
#define SLBM_FORT_SHELL_EXP
 

Functions

SLBM_FORT_SHELL_EXP_IMP void slbm_get_version (char *version, int *size, int *err)
 Retrieve the SLBM Version number. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_error_message (char *errorMessage, int *size)
 Retrieve error message. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_create (int *err)
 Instantiate a SLBM Interface object. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_create_fixed_earth_radius (double *radius, int *err)
 Instantiate a SLBM Interface object with fixed earth radius. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_load_velocity_model (const char *modelFileName, int *nameLength, int *err)
 Load the velocity model into memory from the specified file or directory. This method automatically determines the format of the model. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_save_velocity_model (const char *modelFileName, int *nameLength, int *err)
 Save the velocity model currently in memory to the specified file in format 4. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_save_velocity_model_f (const char *modelFileName, int *nameLength, int *format, int *err)
 Save the velocity model currently in memory to the specified file. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_delete (int *err)
 Deletes the SlbmInterface object instantiated with the call slbm_create(). More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_load_velocity_model_binary (const char *modelDirectory, int *nameLength, int *err)
 Load the velocity model into memory from the specified file or directory. This method automatically determines the format of the model. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_specify_output_directory (const char *directoryName, int *dirNameLength, int *err)
 Deprecated. Use saveVelocityModel() instead. Specify the directory into which the model that is currently in memory should be written the next time that saveVelocityModelBinary() is called. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_save_velocity_model_binary (int *err)
 Deprecated. Use saveVelocityModel() instead. Write the model in format 3 to directory previously specified with a call to specifyOutputDirectory(). More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_create_great_circle (int *phase, double *sourceLat, double *sourceLon, double *sourceDepth, double *receiverLat, double *receiverLon, double *receiverDepth, int *err)
 Instantiate a new GreatCircle object between two locations. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_is_valid (int *err)
 Returns true if the current GreatCirlce object has been instantiated and is ready to be interrogated. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_clear (int *err)
 Delete the current GreatCircle object from memory and clear the pool of stored CrustalProfile objects. The model Grid is not deleted and remains accessible. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_distance (double *dist, int *err)
 Retrieve the source-receiver separation, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_source_distance (double *dist, int *err)
 Retrieve horizontal offset below the source, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_receiver_distance (double *dist, int *err)
 Retrieve horizontal offset below the receiver, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_head_wave_distance (double *dist, int *err)
 Retrieve angular distance traveled by the ray below the headwave interface, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_head_wave_distance_km (double *dist, int *err)
 Retrieve horizontal distance traveled by the ray below the headwave interface, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_travel_time (double *travelTime, int *err)
 Retrieve the total travel time for the GreatCircle, in seconds. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_travel_time_components (double *tTotal, double *tSource, double *tReceiver, double *tHeadwave, double *tGradient, int *err)
 Retrieve the total travel time and the 4 components that contribute to it for the current GreatCircle. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_weights (int nodeId[], double weight[], int *nWeights, int *err)
 Retrieve the weight assigned to each grid node that was touched by the GreatCircle. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_weights_source (int nodeids[], double weights[], int *nWeights, int *err)
 Retrieve the node IDs and the interpolation coefficients for the source CrustalProfile. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_weights_receiver (int nodeids[], double weights[], int *nWeights, int *err)
 Retrieve the node IDs and the interpolation coefficients for the receiver CrustalProfile. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_n_grid_nodes (int *numGridNodes, int *err)
 Retrieve the number of Grid nodes in the Earth model. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_n_head_wave_points (int *nHeadWavePoints, int *err)
 Retrieve the number of LayerProfile objects positioned along the head wave interface. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_grid_data (int *nodeId, double *latitude, double *longitude, double *depth, double *pvelocity, double *svelocity, double *gradient, int *err)
 Retrieve the lat (radians), lon (radians), interface depths (km), P and S wave interval velocities (km/sec) and P and S mantle gradient (1/sec) information associated with a specified node in the velocity grid. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_grid_data (int *nodeId, double *depth, double *pvelocity, double *svelocity, double *gradient, int *err)
 Modify the velocity and gradient information associated with a specified node in the Grid. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_data (int *phase, double *path_increment, double sourceDepth[slbm::NLAYERS], double sourceVelocity[slbm::NLAYERS], double receiverDepth[slbm::NLAYERS], double receiverVelocity[slbm::NLAYERS], double headWaveVelocity[MAX_POINTS], double gradient[MAX_POINTS], int *npoints, int *err)
 Retrieve the data required for input to the travel time calculation, for the current GreatCircle object. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_node_info (int nodes[MAX_POINTS][MAX_NODES], double coefficients[MAX_POINTS][MAX_NODES], int *npoints, int nnodes[MAX_POINTS], int *err)
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_interpolated_point (double *lat, double *lon, int *nodeId, double *coefficients, int *nnodes, double *depth, double *pvelocity, double *svelocity, double *pgradient, double *sgradient, int *err)
 Retrieve interpolated data from the earth model at a single specified latitude, longitude. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_initialize_active_nodes (double *latmin, double *lonmin, double *latmax, double *lonmax, int *err)
 Specify the latitude and longitude range in radians for active nodes. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_n_active_nodes (int *nNodes, int *err)
 Retrieve the number of active nodes in the Grid. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_grid_node_id (int *activeNodeId, int *gridNodeId, int *err)
 Retrieve the grid node ID that corresponds to a specified active node ID. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_id (int *gridNodeId, int *activeNodeId, int *err)
 Retrieve the active node ID that corresponds to a specified grid node ID. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_hit_count (int *nodeId, int *hitCount, int *err)
 Retrieve the number of times that node has been 'touched' by a GreatCircle object. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_neighbors (int *nid, int neighbors[], int *nNeighbors, int *err)
 Retrieve the node IDs of the nodes that surround the specified node. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_neighbor_info (int *nid, int neighbors[], double distance[], double azimuth[], int *nNeighbors, int *err)
 Retrieve the node IDs of the nodes that surround the specified node. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_separation (int *node1, int *node2, double *distance, int *err)
 Retrieve the angular separation of two grid nodes, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_azimuth (int *node1, int *node2, double *azimuth, int *err)
 Retrieve the azimuth from grid node1 to grid node2, radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_traveltime_uncertainty (int *phase, double *distance, double *uncertainty, int *err)
 Calculate a traveltime uncertainty value (seconds) as a function of distance in radians for a supported seismic phase (Pn, Sn, Pg, & Lg). More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_tt_uncertainty (double *tt_uncertainty, int *err)
 Retrieve travel time uncertainty in sec. This function will call either the path-dependent or 1D uncertainty, depending on the model file being used. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_tt_uncertainty_randerr (double *tt_uncertainty, int *err)
 Retrieve travel time uncertainty in sec. This function will call either the path-dependent or 1D uncertainty, depending on the model file being used. This function includes randomError in the computation. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_tt_uncertainty1d (double *tt_uncertainty, int *err)
 Retrieve travel time uncertainty in sec. This function will return the non-path-dependent 1D uncertainty regardless of the model in use. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_zhao_parameters (double *Vm, double *Gm, double *H, double *C, double *Cm, int *udSign, int *err)
 Retrieve some of the parameters that contribute to the calculation of of total travel time using the Zhao algorithm. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_weights (int nodeId[], double weight[], int *nWeights, int *err)
 Retrieve the weight assigned to each active node that was touched by the GreatCircle. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_weights_src (int nodeids[], double weights[], int *nweights, int *err)
 Retrieve the active node IDs and the interpolation coefficients for the source CrustalProfile. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_weights_rcvr (int nodeids[], double weights[], int *nweights, int *err)
 Retrieve the active node IDs and the interpolation coefficients for the receiver CrustalProfile. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_neighbors (int *nid, int neighbors[], int *nNeighbors, int *err)
 Retrieve the node IDs of the nodes that surround the specified node. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_neighbor_info (int *nid, int neighbors[], double distance[], double azimuth[], int *nNeighbors, int *err)
 Retrieve the node IDs of the nodes that surround the specified node. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_data (int *nodeId, double *latitude, double *longitude, double depth[slbm::NLAYERS], double pvelocity[slbm::NLAYERS], double svelocity[slbm::NLAYERS], double gradient[2], int *err)
 Retrieve the lat (radians), lon (radians), interface depths (km), P and S wave interval velocities (km/sec) and P and S mantle gradient (1/sec) information associated with a specified active node in the velocity grid. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_active_node_data (int *nodeId, double depth[slbm::NLAYERS], double pvelocity[slbm::NLAYERS], double svelocity[slbm::NLAYERS], double gradient[2], int *err)
 Modify the velocity and gradient information associated with a specified active node in the Grid. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_ch_max (double *chMax, int *err)
 Set the value of chMax. c is the zhao c parameter and h is the turning depth of the ray below the moho. Zhao method only valid for c*h << 1. When c*h > chMax, then slbm will throw an exception. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_ch_max (double *chMax, int *err)
 Retrieve the current value of chMax. c is the zhao c parameter and h is the turning depth of the ray below the moho. Zhao method only valid for c*h << 1. When c*h > chMax, then slbm will throw an exception. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_average_mantle_velocity (int *type, double *velocity, int *err)
 Retrieve the average P or S wave mantle velocity that is specified in the model input file, in km/sec. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_average_mantle_velocity (int *type, double *velocity, int *err)
 Set the average P or S wave mantle velocity that is recorded in the model input file, in km/sec. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_tess_id (char *tessId, int *size, int *err)
 Retrieve the tessellation ID of the model currently in memory. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_fraction_active (double *fractionActive, int *err)
 Retrieve the fraction of the path length of the current GreatCircle object that is within the currently defined active region. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_max_distance (const double *maxDistance, int *err)
 Set the maximum source-receiver separation for Pn/Sn phase, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_max_distance (double *maxDistance, int *err)
 Retrieve the current value for the maximum source-receiver separation, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_max_depth (const double *maxDepth, int *err)
 Set the maximum source depth for Pn/Sn phase, in km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_max_depth (double *maxDepth, int *err)
 Retrieve the current value for the maximum source depth, in km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_pg_lg_components (double *tTotal, double *tTaup, double *tHeadwave, double *pTaup, double *pHeadwave, double *trTaup, double *trHeadwave, int *err)
 Retrieve information about Pg/Lg travel time calculations. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_dtt_dlat (double *dtt_dlat, int *err)
 Retrieve the derivative of travel time wrt to source latitude, in seconds/radian. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_dtt_dlon (double *dtt_dlon, int *err)
 Retrieve the derivative of travel time wrt to source longitude, in seconds/radian. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_dtt_ddepth (double *dtt_ddepth, int *err)
 Retrieve the derivative of travel time wrt to source depth, in seconds/km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_slowness (double *slowness, int *err)
 Retrieve the horizontal slowness, i.e., the derivative of travel time wrt to receiver-source distance, in seconds/radian. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_slowness_uncertainty (int *phase, double *distance, double *uncertainty, int *err)
 Retrieve uncertainty of the slowness, in seconds/radian for specified phase, distance in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_sh_uncertainty (double *slowness_uncertainty, int *err)
 Retrieve uncertainty of the slowness, in seconds/radian using the phase and distance specified in last call to createGreatCircle(). More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_pierce_point_source (double *lat, double *lon, double *depth, int *err)
 Retrieve the latitude and longitude of the moho pierce point below the source, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_pierce_point_receiver (double *lat, double *lon, double *depth, int *err)
 Retrieve the latitude and longitude of the moho pierce point below the receiver, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_locations (double lat[MAX_POINTS], double lon[MAX_POINTS], double depth[MAX_POINTS], int *npoints, int *err)
 Retrieve the latitudes, longitudes and depths of all the profile positions along the moho. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_points (double *aLat, double *aLon, double *bLat, double *bLon, int *npoints, double latitude[MAX_POINTS], double longitude[MAX_POINTS], int *err)
 Retrieve an array of lat, lon points along a great circle path between two specified points, a and b. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_points_on_centers (double *aLat, double *aLon, double *bLat, double *bLon, int *npoints, double latitude[MAX_POINTS], double longitude[MAX_POINTS], int *err)
 Retrieve an array of lat, lon points along a great circle path between two specified points, a and b. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_dist_az (double *aLat, double *aLon, double *bLat, double *bLon, double *distance, double *azimuth, double *naValue, int *err)
 compute distance and azimuth between two points, A and B (all quantities are in radians). More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_move_point (double *aLat, double *aLon, double *distance, double *azimuth, double *bLat, double *bLon, int *err)
 Find point B that is the specified distance and azimuth from point A. All quantities are in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_del_distance (double *delDistance, int *err)
 Retrieve del_distance in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_del_distance (double *delDistance, int *err)
 Set the value of del_distance, radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_del_depth (double *delDepth, int *err)
 Retrieve del_depth in km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_del_depth (double *delDepth, int *err)
 Set the value of del_depth, in km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_rayparameter (double *rayParameter, int *err)
 Retrieve rayParameter in sec/km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_turningradius (double *turningRadius, int *err)
 Retrieve turning radius in km. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_path_increment (double *pathIncrement, int *err)
 Set the desired spacing of great circle nodes along the head wave interface, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_path_increment (double *pathIncrement, int *err)
 Retrieve the current value of the spacing of great circle nodes along the head wave interface, in radians. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_set_interpolator_type (const char *interpolatorType, int *size, int *err)
 Specify the interpolation type to use, either 'linear' or 'natural_neighbor'. More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_get_interpolator_type (char *interpolatorType, int *size, int *err)
 Retrieve the type of interpolator currently in use; either "LINEAR" or "NATUTAL_NEIGHBOR". More...
 
SLBM_FORT_SHELL_EXP_IMP void slbm_tostring (char *tostring, int *verbosity, int *size, int *err)
 

Detailed Description

FILE slbm_F_shell.h

DESCRIPTION Prototypes for the SLBM Fortran Shell Interface

A Fortran Interface to the SLBM C++ library, providing access to all supported functionality.

The slbm_F_shell Interface (written in C++) manages a connection to the C++ SLBM library and provides Fortran compatible functions. To make this work, 4 conditions must be met:
1) The C++ shared objects libslbm.so and libslbmFshell.so must both be accessible via the LD_LIBRARY_PATH.
2) The Fortran application must be linked to libslbmFshell.so
3) The Fortran application must execute the command: slbm_create();
4) The Fortran application must execute the command: slbm_load_velocity_model_binary();

If all of these conditions are successfully met, then the the Fortran application can use the slbm_F_shell interface through the shared object library libslbmFshell.so to access all the methods described in this document.

Almost all methods in the Fortran interface have a parameter int* err that is set to 0 if the method completed successfully or to a positive error code if the method generated some sort of error. Applications should check the value of err after each call. If err is > 0 then method slbm_get_error_message() can be called to retrieve an errorMessage that provides information about the error. Some of the most important error codes are listed in the followng table.

CodeFunctionMessage
106Grid::reaDataBuffererFromFileCould not open velocity model file.
200GreatCircle_Xg constructorPg/Lg not valid because source or receiver is below the Moho.
201GreatCircle_Xn::computeTravelTimeSource-receiver separation exceeds maximum value.
202GreatCircle_Xn::computeTravelTimeSource depth exceeds maximum value.
203GreatCircle_Xn::computeTravelTimeCrustHorizontal offset below the source plus horizontal offset below the receiver is greater than the source-receiver separation
300GreatCircle_Xn::computeTravelTimeCrustnIterations == 10000
301GreatCircle_Xn::computeTravelTimeCrustc*H is greater than ch_max.
302GreatCircle_Xn::computeTravelTimeMantleCould not converge on stable ray parameter.
303GreatCircle_Xn::computeTravelTimeMantlesearch for minimum H failed.
304GreatCircle_Xn::computeTravelTimeMantlec*H > ch_max.
305GreatCircle_Xn::brentToo many iterations.
400GreatCircle_Xg::computeTravelTimecomputeTravelTimeTaup() and computeTravelTimeHeadwave() both returned NA_VALUE.
401GreatCircle_Xn::computeTravelTimeReceiver depth below Moho is illegal.
402GreatCircle_Xn::computeTravelTimeCrustSource is too close to the receiver.

The libslbmFshell.so maintains a C++ Grid object, which manages interaction with the Earth model. The Earth model is loaded by calling the slbm_load_velocity_model_binary(String) method, described below. This Grid object remains in memory until deleted with the command slbm_delete(), or a different Earth model is loaded with another call to slbm_load_velocity_model_binary(String).

libslbmFshell.so also maintains a single instance of a C++ GreatCircle object which is instantiated with a call to slbm_create_great_circle(). Once instantiated, many slbm_F_shell methods can retrieve information from this GreatCircle object, such as slbm_get_travel_time(), slbm_get_weights(), and more. Once instantiated, the GreatCircle can be interrogated until it is replaced with another GreatCircle by a subsequent call to slbm_create_great_circle(), or is deleted by slbm_clear().

The Grid object stores a vector of map objects which associates the phase and Location of a CrustalProfile object with a pointer to the instance of the CrustalProfile. When slbm_create_great_circle() is called with a latitude, longitude and depth which has been used before, the Grid object will return a pointer to the existing CrustalProfile object, thereby enhancing performance. This vector of maps is cleared when slbm_clear() is called. The implications of all this is that applications that loop over many calls to slbm_create_great_circle() will see a performance improvement if slbm_clear() is not called within the loop. However, for problems where the number of sources and/or receivers is so large that memory becomes an issue, applications could call slbm_clear() within the loop to save memory.

All calculations assume the Earth is defined by a GRS80 ellipsoid. For a description of how points along a great circle are calculated see geovectors.pdf

Definition in file slbm_F_shell.h.

Macro Definition Documentation

◆ SLBM_FORT_SHELL_EXP

#define SLBM_FORT_SHELL_EXP

Definition at line 157 of file slbm_F_shell.h.

◆ SLBM_FORT_SHELL_EXP_IMP

#define SLBM_FORT_SHELL_EXP_IMP

Definition at line 156 of file slbm_F_shell.h.

Function Documentation

◆ slbm_clear()

SLBM_FORT_SHELL_EXP_IMP void slbm_clear ( int *  err)

Delete the current GreatCircle object from memory and clear the pool of stored CrustalProfile objects. The model Grid is not deleted and remains accessible.

Delete the current GreatCircle object from memory and clear the pool of stored CrustalProfile objects. The model Grid is not deleted and remains accessible.

The Grid object owned by SlbmInterface stores a vector of map objects which associates the phase and Location of a CrustalProfile object with a pointer to the instance of the CrustalProfile. When createGreatCircle() is called with a latitude, longitude and depth which has been used before, the Grid object will return a pointer to the existing CrustalProfile object, thereby enhancing performance. This vector of maps is cleared when SlbmInterface::clear() is called. The implications of all this is that applications that loop over many calls to createGreatCircle() will see a performance improvement if clear() is not called within the loop. However, for problems with a huge number of sources and or receivers, if memory becomes an issue, applications could call clear() within the loop to save memory.

Parameters
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_create()

SLBM_FORT_SHELL_EXP_IMP void slbm_create ( int *  err)

Instantiate a SLBM Interface object.

Instantiate a SLBM Interface object.

Parameters
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_create_fixed_earth_radius()

SLBM_FORT_SHELL_EXP_IMP void slbm_create_fixed_earth_radius ( double *  radius,
int *  err 
)

Instantiate a SLBM Interface object with fixed earth radius.

Instantiate a SLBM Interface object fixing the earth radius to the double value passed in argument list.

Parameters
radiusconstant value of earth radius, in km (usually 6371.)
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_create_great_circle()

SLBM_FORT_SHELL_EXP_IMP void slbm_create_great_circle ( int *  phase,
double *  sourceLat,
double *  sourceLon,
double *  sourceDepth,
double *  receiverLat,
double *  receiverLon,
double *  receiverDepth,
int *  err 
)

Instantiate a new GreatCircle object between two locations.

Instantiate a new GreatCircle object between two locations.

Parameters
phasethe phase that this GreatCircle is to support. Recognized phases are Pn, Sn, Pg and Lg.
sourceLatthe geographic latitude of the source in radians.
sourceLonthe longitude of source in radians.
sourceDepththe depth of the source in km.
receiverLatthe geographic latitude of the receiver in radians.
receiverLonthe longitude of the receiver in radians.
receiverDepththe depth of the receiver in km.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_delete()

SLBM_FORT_SHELL_EXP_IMP void slbm_delete ( int *  err)

Deletes the SlbmInterface object instantiated with the call slbm_create().

Deletes the SlbmInterface object instaniated with the call slbm_create().

Parameters
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_data()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_data ( int *  nodeId,
double *  latitude,
double *  longitude,
double  depth[slbm::NLAYERS],
double  pvelocity[slbm::NLAYERS],
double  svelocity[slbm::NLAYERS],
double  gradient[2],
int *  err 
)

Retrieve the lat (radians), lon (radians), interface depths (km), P and S wave interval velocities (km/sec) and P and S mantle gradient (1/sec) information associated with a specified active node in the velocity grid.

Retrieve the interface depth, velocity and gradient information associated with a specified active node in the velocity grid.

Parameters
nodeIdthe active node ID of the grid point in the model (zero based index).
latitudethe latitude of the grid node in radians.
longitudethe longitude of the grid node in radians.
depththe depths of all the model interfaces, in km.
pvelocityan array containing the P velocities of all the intervals at the specified grid node, in km/sec.
svelocityan array containing the S velocities of all the intervals at the specified grid node, in km/sec.
gradienta 2-element array containing the P and S velocity gradients in the mantle, in 1/sec.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_id()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_id ( int *  gridNodeId,
int *  activeNodeId,
int *  err 
)

Retrieve the active node ID that corresponds to a specified grid node ID.

Retrieve the active node ID that corresponds to a specified grid node ID.

Parameters
gridNodeId
activeNodeId
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_neighbor_info()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_neighbor_info ( int *  nid,
int  neighbors[],
double  distance[],
double  azimuth[],
int *  nNeighbors,
int *  err 
)

Retrieve the node IDs of the nodes that surround the specified node.

Retrieve the node IDs of the nodes that surround the specified node.

The caller must supply int array neighbors which is dimensioned large enough to hold the maximum number of neighbors that a node can have, which is 8. The actual number of neighbors is returned in nNeighbors.

Parameters
nid
neighbors
distance
azimuth
nNeighbors
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_neighbors()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_neighbors ( int *  nid,
int  neighbors[],
int *  nNeighbors,
int *  err 
)

Retrieve the node IDs of the nodes that surround the specified node.

Retrieve the node IDs of the nodes that surround the specified node.

Parameters
nidthe id of the node whose neighbors are being requested
neighborsarray of node ids that are the neighbors of nid
nNeighborsthe length of neighbors array
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_weights()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_weights ( int  nodeId[],
double  weight[],
int *  nWeights,
int *  err 
)

Retrieve the weight assigned to each active node that was touched by the GreatCircle.

Retrieve the weight assigned to each active node that Retrieve the weight assigned to each node that was touched by the GreatCircle.

A map which associates an instance of a GridProfile object with a double weight is initialized. Then every LayerProfile on the head wave interface between the source and receiver is visited and the angular distance, d, that the ray traveled in the horizontal segment is retrieved. If d > 0, then the neighboring GridProfile objects that contributed to the interpolated value of the LayerProfile are visited. The product of d * R * C is added to the weight associated with that GridProfile object, where R is the radius of the head wave interface for the LayerProfile object being evaluated, and C is the interpolation coefficient for the GridProfile - LayerProfile pair under consideration. Then, all the GridProfile objects in the map are visited, the grid node IDs extracted into int array nodeId, and the weight extracted into double array weight.

Note: Only grid nodes touched by this GreatCircle are included in the output. Each grid node is included only once, even though more than one LayerProfile object may have contributed some weight to it. The sum of all the weights will equal the horizontal distance traveled by the ray along the head wave interface, from the source pierce point to the receiver pierce point, in km.

Parameters
nodeIdthe active node IDs of all the grid nodes touched by the current GreatCircle.
weightthe weights of all the grid nodes touched by the current GreatCircle. Calling application must dimension this array large enough to handle any possible size.
nWeightsthe number of elements in nodeId and weight. Calling application must dimension this array large enough to handle any possible size.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_weights_rcvr()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_weights_rcvr ( int  nodeids[],
double  weights[],
int *  nweights,
int *  err 
)

Retrieve the active node IDs and the interpolation coefficients for the receiver CrustalProfile.

Retrieve the active node IDs and the interpolation coefficients for the receiver CrustalProfile. There will be anywhere from 1 to 5 of each. The sum of the weights will equal 1.

Parameters
nodeids
weights
nweightssize of nodeids and weights
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_active_node_weights_src()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_active_node_weights_src ( int  nodeids[],
double  weights[],
int *  nweights,
int *  err 
)

Retrieve the active node IDs and the interpolation coefficients for the source CrustalProfile.

Retrieve the active node IDs and the interpolation coefficients for the source CrustalProfile. There will be anywhere from 1 to 5 of each. The sum of the weights will equal 1.

Parameters
nodeids
weights
nweightssize of nodeids and weights
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_average_mantle_velocity()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_average_mantle_velocity ( int *  type,
double *  velocity,
int *  err 
)

Retrieve the average P or S wave mantle velocity that is specified in the model input file, in km/sec.

Retrieve the average P or S wave mantle velocity that is specified in the model input file. This value is used in the calculation of the Zhao c parameter.

Parameters
typespecify either BaseObject::PWAVE or BaseObject::SWAVE.
velocitythe P or S wave velocity is returned in this parameter, in km/sec.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_ch_max()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_ch_max ( double *  chMax,
int *  err 
)

Retrieve the current value of chMax. c is the zhao c parameter and h is the turning depth of the ray below the moho. Zhao method only valid for c*h << 1. When c*h > chMax, then slbm will throw an exception.

Retrieve the current value of chMax. c is the zhao c parameter and h is the turning depth of the ray below the moho. Zhao method only valid for c*h << 1. When c*h > chMax, then slbm will throw an exception. This call retrieves global parameter BaseObject::ch_max

Parameters
chMax
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_del_depth()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_del_depth ( double *  delDepth,
int *  err 
)

Retrieve del_depth in km.

Retrieve current value of del_depth, in km. This is the vertical separation between two points used to compute the derivative of travel time with respect to depth.

Parameters
delDepthvertical separation in km.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_del_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_del_distance ( double *  delDistance,
int *  err 
)

Retrieve del_distance in radians.

Retrieve current value of del_distance, in radians. This is the horizontal separation between two points used to compute horizontal slowness and the derivative of travel time with respect to lat and lon.

Parameters
delDistancehorizontal separation in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_dist_az()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_dist_az ( double *  aLat,
double *  aLon,
double *  bLat,
double *  bLon,
double *  distance,
double *  azimuth,
double *  naValue,
int *  err 
)

compute distance and azimuth between two points, A and B (all quantities are in radians).

compute distance and azimuth between two points, A and B (all quantities are in radians). Computed distance will range between 0 and PI and azimuth will range from -PI to PI. If distance is zero, or if A is located at north or south pole, azimuth will be set to naValue.

Parameters
aLatthe latitude of the first specified point, in radians.
aLonthe longitude of the first specified point, in radians.
bLatthe latitude of the second specified point, in radians.
bLonthe longitude of the second specified point, in radians.
distancefrom point A to point B, in radians.
azimuthfrom point A to point B, in radians.
naValuevalue to return if result is invalid, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_distance ( double *  dist,
int *  err 
)

Retrieve the source-receiver separation, in radians.

Retrieve the source-receiver separation, in radians.

Parameters
distthe source-receiver separation is returned in distance. If the GreatCircle is invalid, distance will equal BaseObject::NA_VALUE.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_dtt_ddepth()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_dtt_ddepth ( double *  dtt_ddepth,
int *  err 
)

Retrieve the derivative of travel time wrt to source depth, in seconds/km.

Retrieve the derivative of travel time wrt to source depth, in seconds/km.

Parameters
dtt_ddepththe derivative of travel time wrt to source depth.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_dtt_dlat()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_dtt_dlat ( double *  dtt_dlat,
int *  err 
)

Retrieve the derivative of travel time wrt to source latitude, in seconds/radian.

Retrieve the derivative of travel time wrt to source latitude, in seconds/radian.

Parameters
dtt_dlatthe derivative of travel time wrt to source latitude.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_dtt_dlon()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_dtt_dlon ( double *  dtt_dlon,
int *  err 
)

Retrieve the derivative of travel time wrt to source longitude, in seconds/radian.

Retrieve the derivative of travel time wrt to source longitude, in seconds/radian.

Parameters
dtt_dlonthe derivative of travel time wrt to source longitude.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_error_message()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_error_message ( char *  errorMessage,
int *  size 
)

Retrieve error message.

Retrieves last error message. This method should be called following any function call that returns a "1" error flag. The string is populated with the text obtained from the SLBMException emessage attribute.

Parameters
errorMessagea char pointer to contain the err message. NOTE: user should allocate storage for approximately 500 chars.
sizenumber of characters allocated for errorMessage.

◆ slbm_get_fraction_active()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_fraction_active ( double *  fractionActive,
int *  err 
)

Retrieve the fraction of the path length of the current GreatCircle object that is within the currently defined active region.

Retrieve the fraction of the path length of the current GreatCircle object that is within the currently defined active region.

Parameters
fractionActive
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_great_circle_data()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_data ( int *  phase,
double *  path_increment,
double  sourceDepth[slbm::NLAYERS],
double  sourceVelocity[slbm::NLAYERS],
double  receiverDepth[slbm::NLAYERS],
double  receiverVelocity[slbm::NLAYERS],
double  headWaveVelocity[MAX_POINTS],
double  gradient[MAX_POINTS],
int *  npoints,
int *  err 
)

Retrieve the data required for input to the travel time calculation, for the current GreatCircle object.

Retrieve the data required for input to the travel time calculation, for the current GreatCircle object.

All the one dimensional arrays can be dimensioned to any size that the developer thinks will be large enough to hold all the values that will result.

The 2D arrays, neighbors and coefficients, must be dimensioned to be exactly MAX_POINTS by nCoefficients. nCoefficients is defined in BaseObject.h and is the number of points that define a cell. For quadrilateral cells, nCoefficients = 4, for triangular cells, nCoefficients = 3. At the time of this release, quadrilateral cells are being used and nCoefficients is 4 (but this may change!). MAX_POINTS is defined in this file, and must be as large or larger than the number of points along the head wave interface. For a great circle path where the source-receiver separation is delta and the angular spacing of the nodes along the head wave interface is path_increment, then MAX_POINTS must be at least delta / path_increment. If MAX_POINTS needs to be increased, it must be changed in this file and slbm_F_shell recompiled (SLBM will not need to be recompiled, however).

Parameters
phasethe phase supported by the current GreatCircle. Will be one of 0=Pn, 1=Sn, 2=Pg, 3=Lg.
path_incrementthe actual horizontal separation of the LayerProfile objects along the head wave interface, in radians. The actual separation will be reduced from the value requested in the call to createGreatCircle() in order that some number of equal sized increments will exactly fit between the source and receiver.
sourceDepththe depths of all the model interfaces below the source, in km.
sourceVelocitythe P or S velocity of each interval below the source, in km/sec.
receiverDepththe depths of all the model interfaces below the receiver, in km.
receiverVelocitythe P or S velocity of each interval below the receiver, in km/sec.
headWaveVelocitythe P or S velocity at the center of each horizontal segment between the source and the receiver, in km/sec. The first horizontal segment starts at the source, the last horizontal segment ends at the receiver, and each one is of size path_increment. The head wave velocities are interpolated at the center of each of these horizontal segments, just below the head wave interface.
gradientthe P or S velocity gradient in the mantle at the center of each horizontal segment of the head wave, in 1/sec. For Pg and Lg, the values will be BaseObject::NA_VALUE.
npointsthe number of horizontal increments sampled along the head wave interface. To discover this number before calling this method call getNHeadWavePoints().
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_great_circle_locations()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_locations ( double  lat[MAX_POINTS],
double  lon[MAX_POINTS],
double  depth[MAX_POINTS],
int *  npoints,
int *  err 
)

Retrieve the latitudes, longitudes and depths of all the profile positions along the moho.

Retrieve the latitudes, longitudes and depths of all the profile positions along the moho. Profile positions are located at the center of each segment of the head wave interface between the source and receiver. The first position is located path_increment/2 radians from the source, the last profile position is located path_increment/2 radians from the receiver, and the others are spaced path_increment radians apart.

Parameters
latthe latitude at the center of each headwave segment, in radians.
lonthe longitude at the center of each headwave segment, in radians.
depththe depth of the headwave interface at the center of each headwave segment, in km.
npointsthe number of horizontal increments sampled along the head wave interface.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_great_circle_node_info()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_node_info ( int  nodes[MAX_POINTS][MAX_NODES],
double  coefficients[MAX_POINTS][MAX_NODES],
int *  npoints,
int  nnodes[MAX_POINTS],
int *  err 
)

◆ slbm_get_great_circle_points()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_points ( double *  aLat,
double *  aLon,
double *  bLat,
double *  bLon,
int *  npoints,
double  latitude[MAX_POINTS],
double  longitude[MAX_POINTS],
int *  err 
)

Retrieve an array of lat, lon points along a great circle path between two specified points, a and b.

Retrieve an array of lat, lon points along a great circle path between two specified points. The great circle path between a and b is divided into npoints-1 equal size cells and the computed points are located at the boundaries of those cells. First point will coincide with point a and last point with point b.

Parameters
aLatthe latitude of the first specified point, in radians.
aLonthe longitude of the first specified point, in radians.
bLatthe latitude of the second specified point, in radians.
bLonthe longitude of the second specified point, in radians.
npointsthe desired number of points along the great circle, in radians.
latitudethe latitudes of the points along the great circle, in radians.
longitudethe longitudes of the points along the great circle, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_great_circle_points_on_centers()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_great_circle_points_on_centers ( double *  aLat,
double *  aLon,
double *  bLat,
double *  bLon,
int *  npoints,
double  latitude[MAX_POINTS],
double  longitude[MAX_POINTS],
int *  err 
)

Retrieve an array of lat, lon points along a great circle path between two specified points, a and b.

Retrieve an array of lat, lon points along a great circle path between two specified points. The great circle path between a and b is divided into npoints equal size cells and the computed points are located at the centers of those cells.

Parameters
aLatthe latitude of the first specified point, in radians.
aLonthe longitude of the first specified point, in radians.
bLatthe latitude of the second specified point, in radians.
bLonthe longitude of the second specified point, in radians.
npointsthe desired number of points along the great circle, in radians.
latitudethe latitudes of the points along the great circle, in radians.
longitudethe longitudes of the points along the great circle, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_grid_data()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_grid_data ( int *  nodeId,
double *  latitude,
double *  longitude,
double *  depth,
double *  pvelocity,
double *  svelocity,
double *  gradient,
int *  err 
)

Retrieve the lat (radians), lon (radians), interface depths (km), P and S wave interval velocities (km/sec) and P and S mantle gradient (1/sec) information associated with a specified node in the velocity grid.

Retrieve the interface depth, velocity and gradient information associated with a specified node in the velocity grid.

Parameters
nodeIdthe node ID of the grid point in the model (zero based index).
latitudethe latitude of the grid node in radians.
longitudethe longitude of the grid node in radians.
depththe depths of all the model interfaces, in km.
pvelocityan array containing the P velocities of all the intervals at the specified grid node, in km/sec.
svelocityan array containing the S velocities of all the intervals at the specified grid node, in km/sec.
gradienta 2-element array containing the P and S velocity gradients in the mantle, in 1/sec.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_grid_node_id()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_grid_node_id ( int *  activeNodeId,
int *  gridNodeId,
int *  err 
)

Retrieve the grid node ID that corresponds to a specified active node ID.

Retrieve the grid node ID that corresponds to a specified active node ID.

Parameters
activeNodeId
gridNodeId
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_head_wave_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_head_wave_distance ( double *  dist,
int *  err 
)

Retrieve angular distance traveled by the ray below the headwave interface, in radians.

Retrieve the angular distance traveled by the ray below the headwave interface, in radians. This is the total distance minus the horizontal offsets below the source and receiver. getSourceDistance() + getReceiverDistance() + getHeadwaveDistance() = getDistance().

Parameters
distthe angular distance traveled by the ray below the headwave interface, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_head_wave_distance_km()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_head_wave_distance_km ( double *  dist,
int *  err 
)

Retrieve horizontal distance traveled by the ray below the headwave interface, in radians.

Retrieve horizontal distance traveled by the ray below the headwave interface, in km. This is the sum of path_increment(i) * R(i) where path_increment(i) is the angular distance traveled by the ray in each angular distance increment along the head wave interface, and R(i) is the radius of the head wave interface in that same horizontal increment.

Parameters
distthe horizontal distance traveled by the ray below the headwave interface, in km.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_interpolated_point()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_interpolated_point ( double *  lat,
double *  lon,
int *  nodeId,
double *  coefficients,
int *  nnodes,
double *  depth,
double *  pvelocity,
double *  svelocity,
double *  pgradient,
double *  sgradient,
int *  err 
)

Retrieve interpolated data from the earth model at a single specified latitude, longitude.

Retrieve interpolated data from the earth model at a single specified latitude, longitude.

Parameters
latthe latitude where information is to be interpolated, in radians.
lonthe longitude where information is to be interpolated, in radians.
nodeIdthe nodeIds of the grid nodes that were involved in the interpolation.
coefficientsthe interpolation coefficients that were applied to the information from the neighboring grid nodes.
nnodesthe number of nodes that contributed to interpolation of data. This is the size of arrays nodeId and coefficients.
depththe depths of the tops of the interfaces in the Earth model, in km. There will be one of these for each layer of the model.
pvelocitythe P velocities of each layer of the model, in km/sec.
svelocitythe S velocities of each layer of the model, in km/sec.
pgradientthe mantle P velocity gradient, in 1/sec.
sgradientthe mantle S velocity gradient, in 1/sec.
Returns
true if successful. If not successful, nodeIds are all -1 and all other returned arrays are populated with BaseObject::NA_VALUE.
Parameters
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_interpolator_type()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_interpolator_type ( char *  interpolatorType,
int *  size,
int *  err 
)

Retrieve the type of interpolator currently in use; either "LINEAR" or "NATUTAL_NEIGHBOR".

Parameters
interpolatorType
sizenumber of characters allocated for interpolatorType
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_max_depth()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_max_depth ( double *  maxDepth,
int *  err 
)

Retrieve the current value for the maximum source depth, in km.

Retrieve the current value for the maximum source depth, in km.

Parameters
maxDepth
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_max_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_max_distance ( double *  maxDistance,
int *  err 
)

Retrieve the current value for the maximum source-receiver separation, in radians.

Retrieve the current value for the maximum source-receiver separation, in radians.

Parameters
maxDistance
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_n_active_nodes()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_n_active_nodes ( int *  nNodes,
int *  err 
)

Retrieve the number of active nodes in the Grid.

Retrieve the number of active nodes in the Grid.

Parameters
nNodes
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_n_grid_nodes()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_n_grid_nodes ( int *  numGridNodes,
int *  err 
)

Retrieve the number of Grid nodes in the Earth model.

Retrieve the number of Grid nodes in the Earth model.

Parameters
numGridNodesthe number of elements in the grid
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_n_head_wave_points()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_n_head_wave_points ( int *  nHeadWavePoints,
int *  err 
)

Retrieve the number of LayerProfile objects positioned along the head wave interface.

Retrieve the number of LayerProfile objects positioned along the head wave interface. It is useful to call this method before calling getGreatCircleData() since the value returned by this method will be the number of elements that will be populated in parameters headWaveVelocity[], neighbors[] and coefficients[].

Parameters
nHeadWavePoints
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_node_azimuth()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_azimuth ( int *  node1,
int *  node2,
double *  azimuth,
int *  err 
)

Retrieve the azimuth from grid node1 to grid node2, radians.

Retrieve the azimuth from grid node1 to grid node2, radians.

Parameters
node1
node2
azimuth
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_node_hit_count()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_hit_count ( int *  nodeId,
int *  hitCount,
int *  err 
)

Retrieve the number of times that node has been 'touched' by a GreatCircle object.

Retrieve the number of times that node has been 'touched' by a GreatCircle object.

Parameters
nodeId
hitCount
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_node_neighbor_info()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_neighbor_info ( int *  nid,
int  neighbors[],
double  distance[],
double  azimuth[],
int *  nNeighbors,
int *  err 
)

Retrieve the node IDs of the nodes that surround the specified node.

Retrieve the node IDs of the nodes that surround the specified node.

The caller must supply int array neighbors which is dimensioned large enough to hold the maximum number of neighbors that a node can have, which is 8. The actual number of neighbors is returned in nNeighbors.

Parameters
nid
neighbors
distance
azimuth
nNeighborslength of neighbors array
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_node_neighbors()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_neighbors ( int *  nid,
int  neighbors[],
int *  nNeighbors,
int *  err 
)

Retrieve the node IDs of the nodes that surround the specified node.

Retrieve the node IDs of the nodes that surround the specified node.

The caller must supply int array neighbors which is dimensioned large enough to hold the maximum number of neighbors that a node can have, which is 8. The actual number of neighbors is returned in nNeighbors.

Parameters
nid
neighbors
nNeighborslength of neighbors array
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_node_separation()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_node_separation ( int *  node1,
int *  node2,
double *  distance,
int *  err 
)

Retrieve the angular separation of two grid nodes, in radians.

Retrieve the angular separation of two grid nodes, in radians.

Parameters
node1
node2
distance
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_path_increment()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_path_increment ( double *  pathIncrement,
int *  err 
)

Retrieve the current value of the spacing of great circle nodes along the head wave interface, in radians.

Retrieve the current value of the spacing of great circle nodes along the head wave interface, in radians. The actual spacing will be reduced from the requested value in order that an integral number of equally spaced LayerProfile objects will exactly span the source-receiver separation. The default value is 0.1 degrees.

Parameters
pathIncrementthe current value of the spacing of great circle nodes along the head wave interface, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_pg_lg_components()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_pg_lg_components ( double *  tTotal,
double *  tTaup,
double *  tHeadwave,
double *  pTaup,
double *  pHeadwave,
double *  trTaup,
double *  trHeadwave,
int *  err 
)

Retrieve information about Pg/Lg travel time calculations.

Retrieve information about Pg/Lg travel time calculations. This method only returns useful information when the phase is Pg or Lg. For Pn and Sn, all information is returned as SLBMGlobals::NA_VALUE.

Parameters
tTotalis the total travel time in seconds. It will be exactly equal to the lesser of tTaup or tHeadwave, except that if tTaup is equal to SLBMGlobals::NA_VALUE, then tTotal will equal tHeadwave.
tTaupis the taup travel time in seconds. If this value is equal to SLBMGlobals::NA_VALUE, it means that the taup calculation failed for some reason (shadow zones, etc.).
tHeadwaveis the headwave travel time in secods
pTaupTauP ray parameter.
pHeadwaveheadwave ray parameter.
trTaupis the radius at which the taup ray turned, in km.
trHeadwaveis the radius at which the headwave ray turned, in km.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_pierce_point_receiver()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_pierce_point_receiver ( double *  lat,
double *  lon,
double *  depth,
int *  err 
)

Retrieve the latitude and longitude of the moho pierce point below the receiver, in radians.

Retrieve the latitude and longitude of the moho pierce point below the receiver, in radians. For Pg, Lg an exception is thrown.

Parameters
latthe latitude of the receiver pierce point, in radians.
lonthe longitude of the receiver pierce point, in radians.
depthmoho depth in km below sea level
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_pierce_point_source()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_pierce_point_source ( double *  lat,
double *  lon,
double *  depth,
int *  err 
)

Retrieve the latitude and longitude of the moho pierce point below the source, in radians.

Retrieve the latitude and longitude of the moho pierce point below the source, in radians. For Pg, Lg and sources in the mantle an exception is thrown.

Parameters
latthe latitude of the source pierce point, in radians.
lonthe longitude of the source pierce point, in radians.
depthmoho depth in km below sea level
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_rayparameter()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_rayparameter ( double *  rayParameter,
int *  err 
)

Retrieve rayParameter in sec/km.

Retrieve the ray parameter in sec/km.

Parameters
rayParameterin sec/km
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_receiver_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_receiver_distance ( double *  dist,
int *  err 
)

Retrieve horizontal offset below the receiver, in radians.

Retrieve horizontal offset below the receiver, in radians. This is the angular distance between the location of the receiver and the receiver pierce point where the ray impinged on the headwave interface.

Parameters
distthe horizontal offset below the receiver, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_sh_uncertainty()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_sh_uncertainty ( double *  slowness_uncertainty,
int *  err 
)

Retrieve uncertainty of the slowness, in seconds/radian using the phase and distance specified in last call to createGreatCircle().

Retrieve uncertainty of the slowness, in seconds/radian using the phase and distance specified in last call to createGreatCircle().

Parameters
slowness_uncertaintyuncertainty of the slowness, in seconds/radian.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_slowness()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_slowness ( double *  slowness,
int *  err 
)

Retrieve the horizontal slowness, i.e., the derivative of travel time wrt to receiver-source distance, in seconds/radian.

Retrieve the horizontal slowness, in seconds/radian.

Parameters
slownessthe derivative of travel time wrt to source latitude.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_slowness_uncertainty()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_slowness_uncertainty ( int *  phase,
double *  distance,
double *  uncertainty,
int *  err 
)

Retrieve uncertainty of the slowness, in seconds/radian for specified phase, distance in radians.

Retrieve uncertainty of the slowness, in seconds/radian for specified phase, distance in radians.

Parameters
phasePn = 0, Sn = 1, Pg = 2, Lg = 3
distancethe angular distance in radians
uncertaintythe model uncertainty (seconds) associated with the angular distance and phase
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_source_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_source_distance ( double *  dist,
int *  err 
)

Retrieve horizontal offset below the source, in radians.

Retrieve horizontal offset below the source, in radians. This is the angular distance between the location of the source and the source pierce point where the ray impinged on the headwave interface.

Parameters
distthe horizontal offset below the source, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_tess_id()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_tess_id ( char *  tessId,
int *  size,
int *  err 
)

Retrieve the tessellation ID of the model currently in memory.

Retrieve the tessellation ID of the model currently in memory. Usually 32 characters long.

Parameters
tessId
sizenumber of characters allocated for tessId
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_travel_time()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_travel_time ( double *  travelTime,
int *  err 
)

Retrieve the total travel time for the GreatCircle, in seconds.

Retrieve the total travel time for the GreatCircle, in seconds.

Parameters
travelTimethe total travel time in seconds is returned in travelTime. If the GreatCircle is invalid, travelTime will equal BaseObject::NA_VALUE.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_travel_time_components()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_travel_time_components ( double *  tTotal,
double *  tSource,
double *  tReceiver,
double *  tHeadwave,
double *  tGradient,
int *  err 
)

Retrieve the total travel time and the 4 components that contribute to it for the current GreatCircle.

Retrieve the total travel time and the 4 components that contribute to it for the current GreatCircle. If the greatCircle is invalid, tTotal and all the components will equal BaseObject::NA_VALUE.

Parameters
tTotalthe total travel time, in seconds.
tSourcethe crustal travel time below the source, in seconds.
tReceiverthe crustal travel time below the receiver, in seconds.
tHeadwavethe head wave travel time, in seconds.
tGradientthe Zhao gradient correction term, in seconds. For GreatCircle objects that support Pg and Lg, this is always 0.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_traveltime_uncertainty()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_traveltime_uncertainty ( int *  phase,
double *  distance,
double *  uncertainty,
int *  err 
)

Calculate a traveltime uncertainty value (seconds) as a function of distance in radians for a supported seismic phase (Pn, Sn, Pg, & Lg).

Calculate a traveltime uncertainty value (seconds) as a function of distance in radians for a supported seismic phase (Pn, Sn, Pg, & Lg).

Parameters
phasePn = 0, Sn = 1, Pg = 2, Lg = 3
distancethe angular distance in radians
uncertaintythe model uncertainty (seconds) associated with the angular distance and phase
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_tt_uncertainty()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_tt_uncertainty ( double *  tt_uncertainty,
int *  err 
)

Retrieve travel time uncertainty in sec. This function will call either the path-dependent or 1D uncertainty, depending on the model file being used.

Retrieve travel time uncertainty in sec. This function will call either the path-dependent or 1D uncertainty, depending on the model file being used.

Parameters
tt_uncertaintyuncertainty of the travel time, in seconds..
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_tt_uncertainty1d()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_tt_uncertainty1d ( double *  tt_uncertainty,
int *  err 
)

Retrieve travel time uncertainty in sec. This function will return the non-path-dependent 1D uncertainty regardless of the model in use.

Retrieve travel time uncertainty in sec using the phase and distance specified in last call to getGreatCircle(). This function will return the non-path-dependent 1D uncertainty regardless of the model in use.

Parameters
tt_uncertaintyuncertainty of the travel time, in seconds..
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_tt_uncertainty_randerr()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_tt_uncertainty_randerr ( double *  tt_uncertainty,
int *  err 
)

Retrieve travel time uncertainty in sec. This function will call either the path-dependent or 1D uncertainty, depending on the model file being used. This function includes randomError in the computation.

Retrieve travel time uncertainty in sec. This function will call either the path-dependent or 1D uncertainty, depending on the model file being used. This function includes randomError in the computation.

Parameters
tt_uncertaintyuncertainty of the travel time, in seconds..
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_turningradius()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_turningradius ( double *  turningRadius,
int *  err 
)

Retrieve turning radius in km.

Retrieve turning radius in km

Parameters
turningRadiusin km.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_version()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_version ( char *  version,
int *  size,
int *  err 
)

Retrieve the SLBM Version number.

Retrieve the SLBM Version number.

Parameters
version
sizenumber of characters allocated for version
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_weights()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_weights ( int  nodeId[],
double  weight[],
int *  nWeights,
int *  err 
)

Retrieve the weight assigned to each grid node that was touched by the GreatCircle.

Retrieve the weight assigned to each grid node that was touched by the GreatCircle.

A map which associates an instance of a GridProfile object with a double weight is initialized. Then every LayerProfile on the head wave interface between the source and receiver is visited and the angular distance, d, that the ray traveled in the horizontal segment is retrieved. If d > 0, then the neighboring GridProfile objects that contributed to the interpolated value of the LayerProfile are visited. The product of d * R * C is added to the weight associated with that GridProfile object, where R is the radius of the head wave interface for the LayerProfile object being evaluated, and C is the interpolation coefficient for the GridProfile - LayerProfile pair under consideration. Then, all the GridProfile objects in the map are visited, the grid node IDs extracted into int array nodeId, and the weight extracted into double array weight.

Note: Only grid nodes touched by this GreatCircle are included in the output. Each grid node is included only once, even though more than one LayerProfile object may have contributed some weight to it. The sum of all the weights will equal the horizontal distance traveled by the ray along the head wave interface, from the source pierce point to the receiver pierce point, in km.

Parameters
nodeIdthe node IDs of all the grid nodes touched by the current GreatCircle.
weightthe weights of all the grid nodes touched by the current GreatCircle. Calling application must dimension this array large enough to handle any possible size.
nWeightsthe number of elements in nodeId and weight. Calling application must dimension this array large enough to handle any possible size.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_weights_receiver()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_weights_receiver ( int  nodeids[],
double  weights[],
int *  nWeights,
int *  err 
)

Retrieve the node IDs and the interpolation coefficients for the receiver CrustalProfile.

Retrieve the node IDs and the interpolation coefficients for the receiver CrustalProfile. There will be BaseObject::nCoefficients of each. The sum of the weights will equal 1.

Parameters
nodeids[]the nodeIds that influenced the interpolated values at the receiver.
weights[]the interpolation coefficients applied to the grid nodes that influenced the interpolated values at the receiver.
nWeightsnumber of weights returned in the weights array
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_weights_source()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_weights_source ( int  nodeids[],
double  weights[],
int *  nWeights,
int *  err 
)

Retrieve the node IDs and the interpolation coefficients for the source CrustalProfile.

Retrieve the node IDs and the interpolation coefficients for the source CrustalProfile. There will be BaseObject::nCoefficients of each. The sum of the weights will equal 1.

Parameters
nodeids[]the nodeIds that influenced the interpolated values at the source.
weights[]the interpolation coefficients applied to the grid nodes that influenced the interpolated values at the source.
nWeightsnumber of weights returned in the weights array
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_get_zhao_parameters()

SLBM_FORT_SHELL_EXP_IMP void slbm_get_zhao_parameters ( double *  Vm,
double *  Gm,
double *  H,
double *  C,
double *  Cm,
int *  udSign,
int *  err 
)

Retrieve some of the parameters that contribute to the calculation of of total travel time using the Zhao algorithm.

Retrieve some of the parameters that contribute to the calculation of of total travel time using the Zhao algorithm. This method only returns meaningful results for phases Pn and Sn. For Pg and Lg, all the parameters of type double are returned with values BaseObject::NA_VALUE and udSign is returned with value of -999.

Parameters
Vmthe velocity at the top of the mantle averaged along the Moho between the source and receiver pierce points.
Gmthe velocity gradient at the top of the mantle averaged along the Moho between the source and receiver pierce points.
Hthe turning depth of the ray relative to the Moho
Ca constant whose product with V0 gives the mantle velocity gradient for a flat Earth. V0 is the velocity of the top of the mantle averaged over the whole model.
Cma constant whose product with Vm gives the mantle velocity gradient for a flat Earth.
udSigna value of 0 indicates the source is in the crust. +1 indicates the ray leaves a mantle source in the downgoing direction. -1 indicates the ray leaves a mantle source in an upgoing direction.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_initialize_active_nodes()

SLBM_FORT_SHELL_EXP_IMP void slbm_initialize_active_nodes ( double *  latmin,
double *  lonmin,
double *  latmax,
double *  lonmax,
int *  err 
)

Specify the latitude and longitude range in radians for active nodes.

Specify the latitude and longitude range in radians for active nodes. Active nodes are defined as follows: for each triangle in the tessellation, if any of the 3 nodes that define the triangle is within the latitude longitude range specified by this method, then all 3 nodes are defined to be active nodes. Lats and lons must be specified in radians.

Parameters
latminminimum latitude in radians
lonminminimum longitude in radians
latmaxmaximum latitude in radians
lonmaxmaximum longitude in radians
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_is_valid()

SLBM_FORT_SHELL_EXP_IMP void slbm_is_valid ( int *  err)

Returns true if the current GreatCirlce object has been instantiated and is ready to be interrogated.

Returns true if the current GreatCirlce object has been instantiated and is ready to be interrogated.

Parameters
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_load_velocity_model()

SLBM_FORT_SHELL_EXP_IMP void slbm_load_velocity_model ( const char *  modelFileName,
int *  nameLength,
int *  err 
)

Load the velocity model into memory from the specified file or directory. This method automatically determines the format of the model.

Load the velocity model into memory from the specified file or directory. This method automatically determines the format of the model and hence is able to load all model formats.

Parameters
modelFileNamethe path to the file or directory that contains the model.
nameLength
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_load_velocity_model_binary()

SLBM_FORT_SHELL_EXP_IMP void slbm_load_velocity_model_binary ( const char *  modelDirectory,
int *  nameLength,
int *  err 
)

Load the velocity model into memory from the specified file or directory. This method automatically determines the format of the model.

Load the velocity model into memory from the specified file or directory. This method automatically determines the format of the model and hence is able to load all model formats.

This method is deprecated in SLBM versions 3 and higher and is provided only for backward compatibility with previous versions of SLBM. It simply calls loadVelocityModel(modelPath).

Parameters
modelDirectory
nameLength
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_move_point()

SLBM_FORT_SHELL_EXP_IMP void slbm_move_point ( double *  aLat,
double *  aLon,
double *  distance,
double *  azimuth,
double *  bLat,
double *  bLon,
int *  err 
)

Find point B that is the specified distance and azimuth from point A. All quantities are in radians.

Find point B that is the specified distance and azimuth from point A. All quantities are in radians.

Parameters
aLatthe latitude of the first specified point, in radians.
aLonthe longitude of the first specified point, in radians.
distancefrom point A to point B, in radians.
azimuthfrom point A to point B, in radians.
bLatthe latitude of the second specified point, in radians.
bLonthe longitude of the second specified point, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_save_velocity_model()

SLBM_FORT_SHELL_EXP_IMP void slbm_save_velocity_model ( const char *  modelFileName,
int *  nameLength,
int *  err 
)

Save the velocity model currently in memory to the specified file in format 4.

Save the velocity model currently in memory to the specified file in format 4.

Parameters
modelFileNamethe full or relative path plus file name of the file to which the earth model is to be written.
nameLength
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_save_velocity_model_binary()

SLBM_FORT_SHELL_EXP_IMP void slbm_save_velocity_model_binary ( int *  err)

Deprecated. Use saveVelocityModel() instead. Write the model in format 3 to directory previously specified with a call to specifyOutputDirectory().

This method is deprecated and is provided only for backward compatibility with previous versions of SLBM. Use method saveVelocityModel() instead.

The model is written in format 3 to the directory previously specified with a call to specifyOutputDirectory().

Parameters
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_save_velocity_model_f()

SLBM_FORT_SHELL_EXP_IMP void slbm_save_velocity_model_f ( const char *  modelFileName,
int *  nameLength,
int *  format,
int *  err 
)

Save the velocity model currently in memory to the specified file.

Save the velocity model currently in memory to the specified file or directory.

The following formats are supported:

  1. SLBM version 1 ascii file. All model information is output to a single file in ascii format. This format was available in SLBM version 2, but never used.

  2. SLBM version 2 directory format. Model information is output to a number of different files and directories, mostly in binary format. This format was used almost exclusively in SLBM version 2.

  3. SLBM version 3 directory format. Model information is output to a number of different files and directories, mostly in binary format. This format is very similar to format 2 with the difference being that the model tessellation and values are stored in GeoTess format instead of the custom SLBM format.

  4. SLBM version 3 single-file format. This is the default+preferred format. All model information is written to a single file. If the modelFileName extension is '.ascii' the file is written in ascii format, otherwise it is written in binary format.

See SLBM_Design.pdf in the main documentation directory for detailed information about model output formats.

Models stored in SLBM version 1 and 2 formats (formats 1 and 2) only support linear interpolation. Models stored in SLBM version 3 formats (formats 3 and 4) support both linear and natural neighbor interpolation.

Parameters
modelFileNamethe full or relative path to the file or directory to which the earth model is to be written. For formats 2 and 3, the directory will be created if it does not exist.
nameLength
formatthe desired format of the output. If omitted, defaults to 4: all model information written to a single file.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_active_node_data()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_active_node_data ( int *  nodeId,
double  depth[slbm::NLAYERS],
double  pvelocity[slbm::NLAYERS],
double  svelocity[slbm::NLAYERS],
double  gradient[2],
int *  err 
)

Modify the velocity and gradient information associated with a specified active node in the Grid.

Modify the velocity and gradient information associated with a specified active node in the Grid.

Parameters
nodeIdthe node number of the grid point in the model. (zero based index).
depthan array containing the depths of the tops of the layers, in km
depthan array containing the depths of the tops of the layers, in km
pvelocityan array containing the P velocities of all the intervals at the specified grid node, in km/sec.
svelocityan array containing the S velocities of all the intervals at the specified grid node, in km/sec.
gradienta 2-element array containing the P and S velocity gradients in the mantle, in 1/sec.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_average_mantle_velocity()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_average_mantle_velocity ( int *  type,
double *  velocity,
int *  err 
)

Set the average P or S wave mantle velocity that is recorded in the model input file, in km/sec.

Set the average P or S wave mantle velocity that is specified in the model input file. This value is used in the calculation of the Zhao c parameter.

Parameters
typespecify either BaseObject::PWAVE or BaseObject::SWAVE.
velocitythe P or S wave velocity that is to be set, in km/sec. This value will be stored in the model file, if the model file is written to file by a call to saveVelocityModel() subsequent to a call to this method.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_ch_max()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_ch_max ( double *  chMax,
int *  err 
)

Set the value of chMax. c is the zhao c parameter and h is the turning depth of the ray below the moho. Zhao method only valid for c*h << 1. When c*h > chMax, then slbm will throw an exception.

Set the value of chMax. c is the zhao c parameter and h is the turning depth of the ray below the moho. Zhao method only valid for c*h << 1. When c*h > chMax, then slbm will throw an exception. This call modifies global parameter BaseObject::ch_max

Parameters
chMax
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_del_depth()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_del_depth ( double *  delDepth,
int *  err 
)

Set the value of del_depth, in km.

Set the value of del_depth, in km. This is the vertical separation between two points used to compute the derivative of travel time with respect to depth.

Parameters
delDepthvertical separation in km.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_del_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_del_distance ( double *  delDistance,
int *  err 
)

Set the value of del_distance, radians.

Set the value of del_distance, in radians. This is the horizontal separation between two points used to compute horizontal slowness and the derivative of travel time with respect to lat and lon.

Parameters
delDistancehorizontal separation in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_grid_data()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_grid_data ( int *  nodeId,
double *  depth,
double *  pvelocity,
double *  svelocity,
double *  gradient,
int *  err 
)

Modify the velocity and gradient information associated with a specified node in the Grid.

Modify the velocity and gradient information associated with a specified node in the Grid.

Parameters
nodeIdthe node number of the grid point in the model. (zero based index).
depthan array containing the depths of the tops of all the intervals at the specified grid node, in km.
pvelocityan array containing the P velocities of all the intervals at the specified grid node, in km/sec.
svelocityan array containing the S velocities of all the intervals at the specified grid node, in km/sec.
gradienta 2-element array containing the P and S velocity gradients in the mantle, in 1/sec.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_interpolator_type()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_interpolator_type ( const char *  interpolatorType,
int *  size,
int *  err 
)

Specify the interpolation type to use, either 'linear' or 'natural_neighbor'.

Specify the interpolation type to use, either "LINEAR" or "NATUTAL_NEIGHBOR". When using the old SLBMGrid objects, LINEAR is the only option allowed.

Parameters
interpolatorType
sizenumber of characters allocated for interpolatorType
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_max_depth()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_max_depth ( const double *  maxDepth,
int *  err 
)

Set the maximum source depth for Pn/Sn phase, in km.

Set the maximum source depth for Pn/Sn phase, in km. Source depths greater than the specified value will result in an exception being thrown in createGreatCircle(). Default value is 9999 km.

Parameters
maxDepth
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_max_distance()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_max_distance ( const double *  maxDistance,
int *  err 
)

Set the maximum source-receiver separation for Pn/Sn phase, in radians.

Set the maximum source-receiver separation for Pn/Sn phase, in radians. Source-receiver separations greater than the specified value will result in an exception being thrown in createGreatCircle(). Default value is PI radians.

Parameters
maxDistance
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_set_path_increment()

SLBM_FORT_SHELL_EXP_IMP void slbm_set_path_increment ( double *  pathIncrement,
int *  err 
)

Set the desired spacing of great circle nodes along the head wave interface, in radians.

Set the desired spacing of great circle nodes along the head wave interface, in radians. The actual spacing will be reduced from the requested value in order that an integral number of equally spaced LayerProfile objects will exactly span the source-receiver separation. Defaults to 0.1 degrees if not specified.

Parameters
pathIncrementthe desired spacing of great circle nodes along the head wave interface, in radians.
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_specify_output_directory()

SLBM_FORT_SHELL_EXP_IMP void slbm_specify_output_directory ( const char *  directoryName,
int *  dirNameLength,
int *  err 
)

Deprecated. Use saveVelocityModel() instead. Specify the directory into which the model that is currently in memory should be written the next time that saveVelocityModelBinary() is called.

This method is deprecated and is provided only for backward compatibility with previous versions of SLBM. Use method saveVelocityModel() instead.

Specify the directory into which the model that is currently in memory should be written the next time that saveVelocityModelBinary() is called. The model will be written in format 3, which is the SLBM version 3 binary directory format (GeoTess).

Parameters
directoryNamethe name of the directory where model files are to be written. The directory will be created if it does not exist.
dirNameLength
err"0" if this call was successful, positive error code if this call generated an error.

◆ slbm_tostring()

SLBM_FORT_SHELL_EXP_IMP void slbm_tostring ( char *  tostring,
int *  verbosity,
int *  size,
int *  err 
)