rod::Rod Struct Reference

#include <rod_structure.h>

Public Member Functions

 Rod (int length, int set_rod_no)
 
 Rod (std::string path, int set_rod_no)
 
Rod set_units ()
 
Rod compute_rest_energy ()
 
Rod do_timestep (RngStream rng[])
 
Rod add_force (float force[4], int node_index)
 
Rod pin_node (bool pin_state, int node_index)
 
Rod load_header (std::string filename)
 
Rod load_contents (std::string filename)
 
Rod write_frame_to_file ()
 
Rod write_array (float *array_ptr, int array_len, float unit_scale_factor)
 
Rod write_mat_params_array (float *array_ptr, int array_len, float stretch_scale_factor, float twist_scale_factor, float length_scale_factor)
 
Rod change_filename (std::string new_filename)
 
Rod equilibrate_rod (RngStream rng[])
 
Rod translate_rod (float *r, float translation_vec[3])
 
Rod rotate_rod (float euler_angles[3])
 
Rod scale_rod (float scale)
 
Rod get_centroid (float *r, float centroid[3])
 
Rod get_min_max (float *r, OUT float min[3], float max[3])
 
Rod get_p (int index, OUT float p[3], bool equil)
 

Data Fields

int length
 
int num_elements
 
int num_rods
 
int rod_no
 
int line_start = 0
 
double rod_version = 999.999
 
bool computed_rest_energy = false
 
float viscosity = 1
 
float timestep = 0.002
 
float kT = 0
 
float perturbation_amount = 0.01
 
float translational_friction
 
float rotational_friction
 
float * equil_r
 
float * equil_m
 
float * current_r
 
float * current_m
 
float * perturbed_x_energy_positive
 
float * perturbed_y_energy_positive
 
float * perturbed_z_energy_positive
 
float * twisted_energy_positive
 
float * perturbed_x_energy_negative
 
float * perturbed_y_energy_negative
 
float * perturbed_z_energy_negative
 
float * twisted_energy_negative
 
float * material_params
 
float * B_matrix
 
float * applied_forces
 
bool * pinned_nodes
 
bool interface_at_start = false
 
bool interface_at_end = false
 
bool restarting = false
 
float bending_response_factor
 
float spring_constant_factor
 
float twist_constant_factor
 
float viscosity_constant_factor
 
std::string rod_filename
 
FILE * file_ptr
 
int frame_no = 0
 
int step_no = 0
 

Constructor & Destructor Documentation

◆ Rod() [1/2]

rod::Rod::Rod ( int  length,
int  set_rod_no 
)

Rod –— Create a rod of a known length. The rod_no is an arbitrary identifier. Note that this creates a rod without initialising the contents of the arrays.

◆ Rod() [2/2]

rod::Rod::Rod ( std::string  path,
int  set_rod_no 
)

Create a rod from a file. This won't do anything other than create the object - it won't even set up any arrays, because it doesn't know how long they'll have to be. After this, you have to call load_header and load_contents, which actually do the dirty work.When we initialize from a file, we don't allocate arrays until we've loaded the file.

Member Function Documentation

◆ add_force()

Rod rod::Rod::add_force ( float  force[4],
int  node_index 
)

Add a constant force that acts upon the rod every timestep. Parameters: force[4], a 4-element array that specifies a 3-D force vector, with the last element being a torque instead. Node_index: the index of the node to apply these forces on. Note: to remove the force, you must call this function again with 0s, or it will continue appyling the force.

Referenced by ffea_test::euler_beam().

◆ change_filename()

Rod rod::Rod::change_filename ( std::string  new_filename)

Close the previous file and create a new file, assigning that to the rod variable *file_ptr. This will also copy the contents of the previous file into this one.

Check if output file exists

Get contents of current file

Create new file

Update member variables (filename, pointer etc)

Referenced by ffea_test::euler_beam(), and World::rod_from_block().

◆ compute_rest_energy()

Rod rod::Rod::compute_rest_energy ( )

◆ do_timestep()

Rod rod::Rod::do_timestep ( RngStream  rng[])

Updates ---—— Do a timestep. This function contains two loops. Both are over the nodes. The first loop populates the contents of the energy arrays, which int we use to work out delta E. The second one uses those energies to compute dynamics and applies those dynamics to the position arrays.

Referenced by ffea_test::connection_propagation(), equilibrate_rod(), ffea_test::euler_beam(), World::run(), and ffea_test::twist_bend_coil().

◆ equilibrate_rod()

Rod rod::Rod::equilibrate_rod ( RngStream  rng[])

Run the simulation for an arbitrary amount of time. If you start a rod exactly in its equilibrium state, chances are it's not going to be equilibrated, which can throw off some tests. It runs for a totally arbitrary 1e-7 seconds and does not save the trajectory from the equilibration.

◆ get_centroid()

Rod rod::Rod::get_centroid ( float *  r,
float  centroid[3] 
)

Get a centroid for the current frame of the rod. Note: you must supply an array (either current_r or equil_r).

Referenced by World::get_system_centroid().

◆ get_min_max()

Rod rod::Rod::get_min_max ( float *  r,
OUT float  min[3],
float  max[3] 
)

◆ get_p()

◆ load_contents()

Rod rod::Rod::load_contents ( std::string  filename)

Load the current state of the rod. This function expects that load_header has already been called. This populates all of the already-initialised arrays containing the state of the rod. Note that it only contains the current state of the rod - the FFEA_rod python class is the only one that loads the rod trajectory.

Make sure this method isn't called before loading header info

Get the line number of the last frame

Check that we got the last frame

Seek back to the start of the file

Find the last frame from the line number we got earlier

Convert each string into a vector<string>

Then convert that into a vector<float>

Check we're not going to overflow and ruin someone's life when we write into the array

Set our rod data arrays to the raw .data() from the vector.

Referenced by ffea_test::euler_beam(), World::rod_from_block(), and ffea_test::twist_bend_coil().

◆ load_header()

Rod rod::Rod::load_header ( std::string  filename)

IO -— Load the header info from a .rodtraj file (that's everything before the —END HEADER— line). Not all the info is read, some of it is for clarity. This populates some rod variables:

  • length - total length of the array (normally 3x the number of nodes)
  • num_elements - number of nodes in the rod
  • num_rods - number of rods in the simulation. Not used right now.
  • line_start - number of the line at which the trajectory begins. This variable is used by load_contents later on, to skip the header.
  • version - which version of the algorithm is this made for? This method will also allocate the memory for all the arrays in the rod. Descriptions of those array are in rod_structure.h. Finally, it sets some default values for global simulation parameters. Eventually, these will be overwritten by parameters from the .ffea file.

This string denotes where the header info ends

Check that we can load our input file

Iterate through file up to rod connections marker.

Use this to prevent rods with bad header info from being initialized

Check that format is valid FFEA_rod

Extract data from lines and deposit it into object

Read in the contents

Set this variable so other methods know we've read the header

Warn the user if there is a file version mismatch

Now that we know the rod length, we can allocate the memory

poiseuille

Referenced by ffea_test::euler_beam(), World::rod_from_block(), and ffea_test::twist_bend_coil().

◆ pin_node()

Rod rod::Rod::pin_node ( bool  pin_state,
int  node_index 
)

◆ rotate_rod()

Rod rod::Rod::rotate_rod ( float  euler_angles[3])

Rotates the rod by the euler angles alpha, beta and gamma (or x, y and z if you prefer. This will update all the node positions AND the material frames will rotate as well. The rotations happen relative to each centroid, so if current_r and equil_r have different centroids, they will be rotated about different points.

Put rod centroid on 0,0,0

Construct rotation matrix from euler angles

Apply rotation matrix

do not try to improve this

Move centroids back

Referenced by World::rod_from_block().

◆ scale_rod()

Rod rod::Rod::scale_rod ( float  scale)

Scale the rod by a float. No return values, it just updates the arrays current_r and equil_r. It doesn't modify m, that'll be normalized away anyway.

Referenced by World::rod_from_block().

◆ set_units()

Rod rod::Rod::set_units ( )

The contents of the rod, by default, are specified in SI units. Although it's possible to do everything in SI, you'll get more precision out of FFEA units. This function will convert all the units into FFEA units. When the file is written, the units are converted back automagically. The units are specified in mesoDimensions.h.

Translate our units into the units specified in FFEA's mesoDimensions header file

And now the rod itself

Return a pointer to the object itself instead of void. Allows for method chaining!

Referenced by ffea_test::euler_beam(), World::rod_from_block(), and ffea_test::twist_bend_coil().

◆ translate_rod()

Rod rod::Rod::translate_rod ( float *  r,
float  translation_vec[3] 
)

Translate every node in the rod by a given translation vector, translation_vec. The parameter float* r is the pointer to any array of node positions, e.g. this->current_r or this->equil_r. No return values, it just updates those arrays

Referenced by World::init(), World::rod_from_block(), and rotate_rod().

◆ write_array()

Rod rod::Rod::write_array ( float *  array_ptr,
int  array_len,
float  unit_scale_factor 
)

Write a single array to a file in the CSV format. Parameters:

  • *array_ptr - the pointer to the array that is to be written.
  • array_len - the length of the array.
  • unit_scale_factor - the unit conversion from the internal FFEA units to SI units.

Referenced by write_frame_to_file().

◆ write_frame_to_file()

Rod rod::Rod::write_frame_to_file ( )

Write the current state of the rod to a file specified by the pointer file_ptr. This will convert from MesoDimensions to SI units, so if your values are already in SI units, they'll be wrong.

Referenced by ffea_test::euler_beam(), World::init(), World::print_trajectory_and_measurement_files(), and ffea_test::twist_bend_coil().

◆ write_mat_params_array()

Rod rod::Rod::write_mat_params_array ( float *  array_ptr,
int  array_len,
float  stretch_scale_factor,
float  twist_scale_factor,
float  length_scale_factor 
)

This function is almost identical to the one above, but it appllies different scale factors for objects in the array,

Referenced by write_frame_to_file().

Field Documentation

◆ applied_forces

float* rod::Rod::applied_forces

Contents of the bending modulus matrix for each node, as a 1-d array. Given as [a_1_1, a_1_2, a_2,1, a_2_2, a_1_1...].

Referenced by add_force(), do_timestep(), and load_header().

◆ B_matrix

float* rod::Rod::B_matrix

◆ bending_response_factor

float rod::Rod::bending_response_factor

If this is true, the rod will skip writing a frame of the trajectory (this is normally done so that the trajectory starts with correct box positioning) Unit conversion factors - the input files are in SI, but internally it uses FFEA's units as determined in dimensions.h

Referenced by set_units(), and write_frame_to_file().

◆ computed_rest_energy

bool rod::Rod::computed_rest_energy = false

If version number unknown, assume latest

◆ current_m

◆ current_r

◆ equil_m

◆ equil_r

float* rod::Rod::equil_r

these will have to be changed when I end up implementing per-element radius, to be computed on-the-fly instead most likely Each set of rod data is stored in a single, c-style array that goes {x,y,z,x,y,z...}

Referenced by ffea_test::connection_propagation(), do_timestep(), rod::Rod_blob_interface::get_attachment_material_axis(), rod::Rod_blob_interface::get_node_energy(), get_p(), load_contents(), load_header(), rod::Rod_blob_interface::position_rod_from_blob(), World::rod_from_block(), rotate_rod(), scale_rod(), set_units(), and write_frame_to_file().

◆ file_ptr

◆ frame_no

int rod::Rod::frame_no = 0

Referenced by do_timestep(), and write_frame_to_file().

◆ interface_at_end

bool rod::Rod::interface_at_end = false

if this is true, the positioning of the start node is being handled by a rod-blob interface, so its energies will be ignored.

Referenced by rod::Rod_blob_interface::Rod_blob_interface().

◆ interface_at_start

bool rod::Rod::interface_at_start = false

This array is the length of the number of nodes in the rod, and it contains a boolean stating whether that node is pinned (true) or not (false).

Referenced by rod::Rod_blob_interface::Rod_blob_interface().

◆ kT

◆ length

◆ line_start

int rod::Rod::line_start = 0

Each rod will be created with a unique ID

Referenced by load_contents(), and load_header().

◆ material_params

◆ num_elements

◆ num_rods

int rod::Rod::num_rods

The number of elements in the rod

Referenced by load_header().

◆ perturbation_amount

float rod::Rod::perturbation_amount = 0.01

◆ perturbed_x_energy_negative

float* rod::Rod::perturbed_x_energy_negative

◆ perturbed_x_energy_positive

float* rod::Rod::perturbed_x_energy_positive

Current configuration of the material frame.

Referenced by do_timestep(), load_contents(), load_header(), set_units(), and write_frame_to_file().

◆ perturbed_y_energy_negative

float* rod::Rod::perturbed_y_energy_negative

Stretch, bend, twist, stretch, bend, twist...

Referenced by do_timestep(), load_contents(), load_header(), set_units(), and write_frame_to_file().

◆ perturbed_y_energy_positive

float* rod::Rod::perturbed_y_energy_positive

Energies associated with the perturbations we do in order to get dE. Given as [stretch, bend, twist, stretch, bend, twist...]

Referenced by do_timestep(), load_contents(), load_header(), set_units(), and write_frame_to_file().

◆ perturbed_z_energy_negative

float* rod::Rod::perturbed_z_energy_negative

◆ perturbed_z_energy_positive

float* rod::Rod::perturbed_z_energy_positive

◆ pinned_nodes

bool* rod::Rod::pinned_nodes

Another [x,y,z,x,y,z...] array, this one containing the force vectors acting on each node in the rod.

Referenced by ffea_test::connection_propagation(), do_timestep(), ffea_test::euler_beam(), load_header(), and ffea_test::twist_bend_coil().

◆ restarting

bool rod::Rod::restarting = false

if this is true, the positioning of the end node is being handled by a rod-blob interface, so its energies will be ignored.

Referenced by World::print_trajectory_and_measurement_files(), and World::rod_from_block().

◆ rod_filename

std::string rod::Rod::rod_filename

Referenced by change_filename(), and load_header().

◆ rod_no

int rod::Rod::rod_no

When the system is set up in ffeatools, this will be set

Referenced by Rod(), World::rod_from_block(), and write_frame_to_file().

◆ rod_version

double rod::Rod::rod_version = 999.999

Keeps track of whether the header file has been read

Referenced by load_contents(), and load_header().

◆ rotational_friction

float rod::Rod::rotational_friction

Referenced by do_timestep().

◆ spring_constant_factor

float rod::Rod::spring_constant_factor

Referenced by set_units(), and write_frame_to_file().

◆ step_no

int rod::Rod::step_no = 0

Referenced by do_timestep().

◆ timestep

◆ translational_friction

float rod::Rod::translational_friction

Amount by which nodes are perturbed during numerical differentiation. May want to override with a local value depending on the scale of the simulation.

Referenced by do_timestep().

◆ twist_constant_factor

float rod::Rod::twist_constant_factor

Referenced by set_units(), and write_frame_to_file().

◆ twisted_energy_negative

float* rod::Rod::twisted_energy_negative

◆ twisted_energy_positive

float* rod::Rod::twisted_energy_positive

◆ viscosity

float rod::Rod::viscosity = 1

Global simulation parameters - eventually we may read this stuff in from the .ffea file

Referenced by do_timestep(), load_header(), World::read_and_build_system(), and ffea_test::twist_bend_coil().

◆ viscosity_constant_factor

float rod::Rod::viscosity_constant_factor

Referenced by load_header().


The documentation for this struct was generated from the following files: