ASPECT
Public Member Functions | Private Member Functions | Private Attributes | List of all members
aspect::BoundaryVelocity::internal::GPlatesLookup< dim > Class Template Reference

Public Member Functions

 GPlatesLookup (const Tensor< 1, 2 > &surface_point_one, const Tensor< 1, 2 > &surface_point_two)
 
std::string screen_output (const Tensor< 1, 2 > &surface_point_one, const Tensor< 1, 2 > &surface_point_two) const
 
void load_file (const std::string &filename, const MPI_Comm comm)
 
Tensor< 1, dim > surface_velocity (const Point< dim > &position) const
 

Private Member Functions

std::array< double, 3 > angles_from_matrix (const Tensor< 2, 3 > &rotation_matrix) const
 
double rotation_axis_from_matrix (Tensor< 1, 3 > &rotation_axis, const Tensor< 2, 3 > &rotation_matrix) const
 
template<int in, int out>
Tensor< 1, out > convert_tensor (const Tensor< 1, in > &old_tensor) const
 
Tensor< 1, 3 > cartesian_surface_coordinates (const Tensor< 1, 3 > &sposition) const
 
Tensor< 1, dim > cartesian_velocity_at_surface_point (const std::array< double, 3 > &spherical_point) const
 
Tensor< 1, 3 > sphere_to_cart_velocity (const Tensor< 1, 2 > &s_velocities, const std::array< double, 3 > &s_position) const
 
bool gplates_1_4_or_higher (const boost::property_tree::ptree &pt) const
 

Private Attributes

std::array< std::unique_ptr< Functions::InterpolatedUniformGridData< 2 > >, 2 > velocities
 
double delta_phi
 
double delta_theta
 
Tensor< 2, 3 > rotation_matrix
 

Detailed Description

template<int dim>
class aspect::BoundaryVelocity::internal::GPlatesLookup< dim >

GPlatesLookup handles all kinds of tasks around looking up a certain velocity boundary condition from a gplates .gpml file.

Definition at line 45 of file gplates.h.

Constructor & Destructor Documentation

§ GPlatesLookup()

template<int dim>
aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::GPlatesLookup ( const Tensor< 1, 2 > &  surface_point_one,
const Tensor< 1, 2 > &  surface_point_two 
)

Initialize all members and calculates any necessary rotation parameters for a 2D model.

Member Function Documentation

§ screen_output()

template<int dim>
std::string aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::screen_output ( const Tensor< 1, 2 > &  surface_point_one,
const Tensor< 1, 2 > &  surface_point_two 
) const

Outputs the GPlates module information at model start.

§ load_file()

template<int dim>
void aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::load_file ( const std::string &  filename,
const MPI_Comm  comm 
)

Loads a gplates .gpml velocity file. Throws an exception if the file does not exist.

§ surface_velocity()

template<int dim>
Tensor<1,dim> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::surface_velocity ( const Point< dim > &  position) const

Returns the computed surface velocity in cartesian coordinates. Takes as input the position. Actual velocity interpolation is performed in spherical coordinates.

Parameters
positionThe current position to compute velocity

§ angles_from_matrix()

template<int dim>
std::array<double,3> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::angles_from_matrix ( const Tensor< 2, 3 > &  rotation_matrix) const
private

A function that returns the corresponding paraview angles for a rotation described by a rotation matrix. These differ from the usually used euler angles by assuming a rotation around the coordinate axes in the order y-x-z (instead of the often used z-x-z)

§ rotation_axis_from_matrix()

template<int dim>
double aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::rotation_axis_from_matrix ( Tensor< 1, 3 > &  rotation_axis,
const Tensor< 2, 3 > &  rotation_matrix 
) const
private

A function that returns the corresponding rotation axis/angle for a rotation described by a rotation matrix.

§ convert_tensor()

template<int dim>
template<int in, int out>
Tensor<1,out> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::convert_tensor ( const Tensor< 1, in > &  old_tensor) const
private

Convert a tensor of rank 1 and dimension in to rank 1 and dimension out. If \(out < in\) the last elements will be discarded, if \(out > in\) zeroes will be appended to fill the tensor.

§ cartesian_surface_coordinates()

template<int dim>
Tensor<1,3> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::cartesian_surface_coordinates ( const Tensor< 1, 3 > &  sposition) const
private

Return the cartesian coordinates of a spherical surface position defined by theta (polar angle, not geographical latitude) and phi.

§ cartesian_velocity_at_surface_point()

template<int dim>
Tensor<1,dim> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::cartesian_velocity_at_surface_point ( const std::array< double, 3 > &  spherical_point) const
private

This function looks up the north- and east-velocities at a given position and converts them to cartesian velocities.

§ sphere_to_cart_velocity()

template<int dim>
Tensor<1,3> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::sphere_to_cart_velocity ( const Tensor< 1, 2 > &  s_velocities,
const std::array< double, 3 > &  s_position 
) const
private

Returns cartesian velocities calculated from surface velocities and position in spherical coordinates

Parameters
s_velocitiesSurface velocities in spherical coordinates (theta, phi)
s_positionPosition in spherical coordinates (radius,phi,theta)

§ gplates_1_4_or_higher()

template<int dim>
bool aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::gplates_1_4_or_higher ( const boost::property_tree::ptree &  pt) const
private

Check whether the gpml file was created by GPlates1.4 or later. We need to know this, because the mesh has changed its longitude origin from 0 to -180 degrees and we need to correct for this.

Member Data Documentation

§ velocities

template<int dim>
std::array<std::unique_ptr<Functions::InterpolatedUniformGridData<2> >, 2> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::velocities
private

Interpolation functions to access the velocities.

Definition at line 83 of file gplates.h.

§ delta_phi

template<int dim>
double aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::delta_phi
private

Distances between adjacent point in the Lat/Long grid

Definition at line 88 of file gplates.h.

§ delta_theta

template<int dim>
double aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::delta_theta
private

Definition at line 88 of file gplates.h.

§ rotation_matrix

template<int dim>
Tensor<2,3> aspect::BoundaryVelocity::internal::GPlatesLookup< dim >::rotation_matrix
private

The matrix, which describes the rotation by which a 2D model needs to be transformed to a plane that contains the origin and the two user prescribed points. Is not necessary and therefore not used for 3D models.

Definition at line 96 of file gplates.h.


The documentation for this class was generated from the following file: