Namespaces | |
detail | |
Classes | |
struct | GravmagInterpolatorParams |
Structure containing gravity and magnetic interpolation parameters. More... | |
Functions | |
template<> | |
ContactPointCache | generateCache< ForwardModel::CONTACTPOINT > (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const ContactPointSpec &spec) |
Generate a cache object for a contact point forward model. More... | |
template<> | |
ContactPointResults | forwardModel< ForwardModel::CONTACTPOINT > (const ContactPointSpec &spec, const ContactPointCache &cache, const WorldParams &world) |
Run a contact point forward model. More... | |
template<ForwardModel f> | |
Types< f >::Cache | generateCache (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const typename Types< f >::Spec &spec) |
Generate a cache object for a specific forward model. Cache objects contain repeatly used information that only needs to be computed once by the forward model. More... | |
template<ForwardModel f> | |
Types< f >::Results | forwardModel (const typename Types< f >::Spec &spec, const typename Types< f >::Cache &cache, const WorldParams &world) |
Run a particular forward model. More... | |
GlobalCache | generateGlobalCache (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const GlobalSpec &spec, const std::set< ForwardModel > &enabled) |
Generate a global cache object for all forward models. More... | |
GlobalResults | forwardModelAll (const GlobalSpec &spec, const GlobalCache &cache, const GlobalParams ¶ms, const std::set< ForwardModel > &enabled) |
Run all the forward models. More... | |
template<> | |
GravCache | generateCache< ForwardModel::GRAVITY > (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const GravSpec &gravSpec) |
Generate a cache object for a gravity forward model. More... | |
template<> | |
GravResults | forwardModel< ForwardModel::GRAVITY > (const GravSpec &spec, const GravCache &cache, const WorldParams &world) |
Run a gravity forward model. More... | |
Eigen::MatrixXd | gravSens (const Eigen::VectorXd &xEdges, const Eigen::VectorXd &yEdges, const Eigen::VectorXd &zEdges, const Eigen::MatrixXd &locations) |
Compute the gravity sensitivity matrix. More... | |
GravmagInterpolatorParams | makeInterpParams (const VoxelSpec &voxelisation, const Eigen::MatrixXd &locations, const WorldSpec &worldSpec) |
Create a GravmagInterpolatorParams object. More... | |
Eigen::VectorXd | computeField (const Eigen::MatrixXd &sens, const Eigen::MatrixXi sensorIndices, const Eigen::MatrixXd sensorWeights, const Eigen::VectorXd &properties) |
Computes the field values for either gravity or magnetic. More... | |
template<> | |
MagCache | generateCache< ForwardModel::MAGNETICS > (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const MagSpec &magSpec) |
Generate a cache object for a magnetic forward model. More... | |
template<> | |
MagResults | forwardModel< ForwardModel::MAGNETICS > (const MagSpec &spec, const MagCache &cache, const WorldParams &world) |
Run a magnetic forward model. More... | |
Eigen::MatrixXd | magSens (const Eigen::VectorXd &xEdges, const Eigen::VectorXd &yEdges, const Eigen::VectorXd &zEdges, const Eigen::MatrixXd &locations, const double &bX, const double &bY, const double &bZ) |
Compute the magnetic sensitivity matrix. More... | |
Eigen::MatrixXd | magSens (const Eigen::VectorXd &xEdges, const Eigen::VectorXd &yEdges, const Eigen::VectorXd &zEdges, const Eigen::MatrixXd &locations, const Eigen::VectorXd &bX, const Eigen::VectorXd &bY, const Eigen::VectorXd &bZ) |
Compute the magnetic sensitivity matrix. More... | |
Eigen::MatrixX4d | phaseTensor1d (const Eigen::MatrixX4cd &Z) |
Calculate the phase tensor given the impedance matrix. More... | |
Eigen::VectorXd | alpha1d (const Eigen::MatrixX4d &tensor) |
Computes the alpha value of the phase tensor, which is a measure of its rotations from the major axis. More... | |
Eigen::VectorXd | beta1d (const Eigen::MatrixX4d &tensor) |
Computes the beta value of the phase tensor, which is a measure of its rotations from the z-axis. More... | |
Eigen::MatrixX4cd | impedenceAniso1d (const Eigen::VectorXd &thicknesses, const Eigen::VectorXd &freqs, const Eigen::VectorXd &resx, const Eigen::VectorXd &resy, const Eigen::VectorXd &phases) |
Calculate the impedence matrix for anisotropic 1D MT. More... | |
template<> | |
MtAnisoCache | generateCache< ForwardModel::MTANISO > (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const MtAnisoSpec &mtSpec) |
Generate a cache object for a MT forward model. More... | |
template<> | |
MtAnisoResults | forwardModel< ForwardModel::MTANISO > (const MtAnisoSpec &spec, const MtAnisoCache &cache, const WorldParams &world) |
Run a MT forward model. More... | |
template<> | |
Seismic1dCache | generateCache< ForwardModel::SEISMIC1D > (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const Seismic1dSpec &spec) |
Generate a cache object for a seismic forward model. More... | |
template<> | |
Seismic1dResults | forwardModel< ForwardModel::SEISMIC1D > (const Seismic1dSpec &spec, const Seismic1dCache &cache, const WorldParams &world) |
Run a seismic forward model. More... | |
boost::multi_array< double, 3 > | fillAndPad (const Eigen::VectorXd &values, const ThermalSpec &spec) |
boost::multi_array< double, 3 > | eastWest (const boost::multi_array< double, 3 > &pad) |
boost::multi_array< double, 3 > | cells (const boost::multi_array< double, 3 > &pad, uint axis) |
boost::multi_array< double, 3 > | isocells (const boost::multi_array< double, 3 > &pad) |
boost::multi_array< double, 3 > | temp (const boost::multi_array< double, 3 > &eastwest, const boost::multi_array< double, 3 > &northsouth, const boost::multi_array< double, 3 > &updown, const boost::multi_array< double, 3 > &sCells, const Eigen::MatrixXd &tempZ0, const Eigen::MatrixXd &zLowBound, bool zLowBoundIsHeatFlow, double xSize, double ySize, double zSize) |
uint | clip (int in, uint maxX) |
Eigen::VectorXd | evalAtLocations (const boost::multi_array< double, 3 > &tempVox, const Eigen::MatrixXd &locations, const std::pair< double, double > &xSize, const std::pair< double, double > &ySize, const std::pair< double, double > &zSize) |
template<> | |
ThermalCache | generateCache< ForwardModel::THERMAL > (const std::vector< world::InterpolatorSpec > &boundaryInterpolation, const WorldSpec &worldSpec, const ThermalSpec &thermSpec) |
template<> | |
ThermalResults | forwardModel< ForwardModel::THERMAL > (const ThermalSpec &spec, const ThermalCache &cache, const WorldParams &world) |
Namespace for all the forward models.
Eigen::VectorXd obsidian::fwd::alpha1d | ( | const Eigen::MatrixX4d & | tensor | ) |
Computes the alpha value of the phase tensor, which is a measure of its rotations from the major axis.
tensor | The phase tensor computed by phaseTensor1d() |
Eigen::VectorXd obsidian::fwd::beta1d | ( | const Eigen::MatrixX4d & | tensor | ) |
Computes the beta value of the phase tensor, which is a measure of its rotations from the z-axis.
tensor | The phase tensor computed by phaseTensor1d() |
Eigen::VectorXd obsidian::fwd::computeField | ( | const Eigen::MatrixXd & | sens, |
const Eigen::MatrixXi | sensorIndices, | ||
const Eigen::MatrixXd | sensorWeights, | ||
const Eigen::VectorXd & | properties | ||
) |
Computes the field values for either gravity or magnetic.
sens | The gravity or magnetic sensitivity matrix. |
sensorIndices | The indices of each of the sensors. |
sensorWeights | The weights of each of the sensors. |
properties | The rock property for the sensor type. |
Types<f>::Results obsidian::fwd::forwardModel | ( | const typename Types< f >::Spec & | spec, |
const typename Types< f >::Cache & | cache, | ||
const WorldParams & | world | ||
) |
Run a particular forward model.
spec | The forward model specification. |
cache | The forward model cache generated by generateCache(). |
world | The world model parameters. |
ContactPointResults obsidian::fwd::forwardModel< ForwardModel::CONTACTPOINT > | ( | const ContactPointSpec & | spec, |
const ContactPointCache & | cache, | ||
const WorldParams & | world | ||
) |
Run a contact point forward model.
spec | The forward model specification. |
cache | The forward model cache generated by generateCache(). |
world | The world model parameters. |
GravResults obsidian::fwd::forwardModel< ForwardModel::GRAVITY > | ( | const GravSpec & | spec, |
const GravCache & | cache, | ||
const WorldParams & | world | ||
) |
Run a gravity forward model.
spec | The forward model specification. |
cache | The forward model cache generated by generateCache(). |
world | The world model parameters. |
MagResults obsidian::fwd::forwardModel< ForwardModel::MAGNETICS > | ( | const MagSpec & | spec, |
const MagCache & | cache, | ||
const WorldParams & | world | ||
) |
Run a magnetic forward model.
spec | The forward model specification. |
cache | The forward model cache generated by generateCache(). |
world | The world model parameters. |
MtAnisoResults obsidian::fwd::forwardModel< ForwardModel::MTANISO > | ( | const MtAnisoSpec & | spec, |
const MtAnisoCache & | cache, | ||
const WorldParams & | world | ||
) |
Run a MT forward model.
spec | The forward model specification. |
cache | The forward model cache generated by generateCache(). |
world | The world model parameters. |
Seismic1dResults obsidian::fwd::forwardModel< ForwardModel::SEISMIC1D > | ( | const Seismic1dSpec & | spec, |
const Seismic1dCache & | cache, | ||
const WorldParams & | world | ||
) |
Run a seismic forward model.
spec | The forward model specification. |
cache | The forward model cache generated by generateCache(). |
world | The world model parameters. |
GlobalResults obsidian::fwd::forwardModelAll | ( | const GlobalSpec & | spec, |
const GlobalCache & | cache, | ||
const GlobalParams & | params, | ||
const std::set< ForwardModel > & | enabled | ||
) |
Run all the forward models.
spec | The global world specifications. |
cache | The global forward model cache generated by generateCacheGlobal(). |
params | The global world parameters. |
Types<f>::Cache obsidian::fwd::generateCache | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const WorldSpec & | worldSpec, | ||
const typename Types< f >::Spec & | spec | ||
) |
Generate a cache object for a specific forward model. Cache objects contain repeatly used information that only needs to be computed once by the forward model.
boundaryInterpolation | The world model interpolation parameters. |
worldSpec | The world model specification. |
spec | The forward model specification. |
ContactPointCache obsidian::fwd::generateCache< ForwardModel::CONTACTPOINT > | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const WorldSpec & | worldSpec, | ||
const ContactPointSpec & | spec | ||
) |
Generate a cache object for a contact point forward model.
boundaryInterpolation | The world model interpolation parameters. |
worldSpec | The world model specification. |
spec | The forward model specification. |
GravCache obsidian::fwd::generateCache< ForwardModel::GRAVITY > | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const WorldSpec & | worldSpec, | ||
const GravSpec & | gravSpec | ||
) |
Generate a cache object for a gravity forward model.
boundaryInterpolation | The world model interpolation parameters. |
worldSpec | The world model specification. |
gravSpec | The forward model specification. |
MagCache obsidian::fwd::generateCache< ForwardModel::MAGNETICS > | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const WorldSpec & | worldSpec, | ||
const MagSpec & | magSpec | ||
) |
Generate a cache object for a magnetic forward model.
boundaryInterpolation | The world model interpolation parameters. |
worldSpec | The world model specification. |
magSpec | The forward model specification. |
MtAnisoCache obsidian::fwd::generateCache< ForwardModel::MTANISO > | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const WorldSpec & | worldSpec, | ||
const MtAnisoSpec & | mtSpec | ||
) |
Generate a cache object for a MT forward model.
boundaryInterpolation | The world model interpolation parameters. |
worldSpec | The world model specification. |
mtSpec | The forward model specification. |
Seismic1dCache obsidian::fwd::generateCache< ForwardModel::SEISMIC1D > | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const WorldSpec & | worldSpec, | ||
const Seismic1dSpec & | spec | ||
) |
Generate a cache object for a seismic forward model.
boundaryInterpolation | The world model interpolation parameters. |
worldSpec | The world model specification. |
spec | The forward model specification. |
GlobalCache obsidian::fwd::generateGlobalCache | ( | const std::vector< world::InterpolatorSpec > & | boundaryInterpolation, |
const GlobalSpec & | spec, | ||
const std::set< ForwardModel > & | enabled | ||
) |
Generate a global cache object for all forward models.
boundaryInterpolation | The world model interpolation parameters. |
spec | The global world specifications. |
Eigen::MatrixXd obsidian::fwd::gravSens | ( | const Eigen::VectorXd & | xEdges, |
const Eigen::VectorXd & | yEdges, | ||
const Eigen::VectorXd & | zEdges, | ||
const Eigen::MatrixXd & | locations | ||
) |
Compute the gravity sensitivity matrix.
xEdges | The x coordinates of the mesh grid. |
yEdges | The y coordinates of the mesh grid. |
zEdges | The z coordinates of the mesh grid. |
locations | A Nx3 matrix containing the coordinates of the sensor locations. |
Eigen::MatrixX4cd obsidian::fwd::impedenceAniso1d | ( | const Eigen::VectorXd & | thicknesses, |
const Eigen::VectorXd & | freqs, | ||
const Eigen::VectorXd & | resx, | ||
const Eigen::VectorXd & | resy, | ||
const Eigen::VectorXd & | phases | ||
) |
Calculate the impedence matrix for anisotropic 1D MT.
thicknesses | A vector containing the thicknesses for each layer. |
freqs | A vector containing the MT frequencies. |
resx | A vector containing the resistivities in x direction (not logarithmic). |
resy | A vector containing the resistivities in y direction (not logarithmic). |
phases | A vector containing the phase of each layer. |
Eigen::MatrixXd obsidian::fwd::magSens | ( | const Eigen::VectorXd & | xEdges, |
const Eigen::VectorXd & | yEdges, | ||
const Eigen::VectorXd & | zEdges, | ||
const Eigen::MatrixXd & | locations, | ||
const double & | bX, | ||
const double & | bY, | ||
const double & | bZ | ||
) |
Compute the magnetic sensitivity matrix.
xEdges | The x coordinates of the mesh grid. |
yEdges | The y coordinates of the mesh grid. |
zEdges | The z coordinates of the mesh grid. |
locations | A Nx3 matrix containing the coordinates of the sensor locations. |
bX,bY,bZ | The field values along each axis. Assuming the field values are constant in each axis. |
Eigen::MatrixXd obsidian::fwd::magSens | ( | const Eigen::VectorXd & | xEdges, |
const Eigen::VectorXd & | yEdges, | ||
const Eigen::VectorXd & | zEdges, | ||
const Eigen::MatrixXd & | locations, | ||
const Eigen::VectorXd & | bX, | ||
const Eigen::VectorXd & | bY, | ||
const Eigen::VectorXd & | bZ | ||
) |
Compute the magnetic sensitivity matrix.
xEdges | The x coordinates of the mesh grid. |
yEdges | The y coordinates of the mesh grid. |
zEdges | The z coordinates of the mesh grid. |
locations | A Nx3 matrix containing the coordinates of the sensor locations. |
bX,bY,bZ | The field values at each of the mesh grid points. |
GravmagInterpolatorParams obsidian::fwd::makeInterpParams | ( | const VoxelSpec & | voxelisation, |
const Eigen::MatrixXd & | locations, | ||
const WorldSpec & | worldSpec | ||
) |
Create a GravmagInterpolatorParams object.
voxelisation | The voxelised world. |
locations | The sensor locations. |
worldSpec | The world specifications. |
Eigen::MatrixX4d obsidian::fwd::phaseTensor1d | ( | const Eigen::MatrixX4cd & | Z | ) |
Calculate the phase tensor given the impedance matrix.
Z | The impedance matrix. |