![]() |
Stellarator-Tools
|
Base class containing all the data needed to reconstruct a model. More...
Public Member Functions | |
FINAL | reconstruction_destruct (this) |
Deconstruct a reconstruction_class object. | |
procedure | get_k_use (this, num_sv) |
Get the number of singular values to use. | |
procedure | get_exp_dg2 (this, delta_a) |
Get the expected dg^2 size. | |
procedure | get_exp_g2 (this, delta_a) |
Get the expected g^2 size. | |
procedure | get_g2 (this) |
Get the current g^2. | |
procedure | get_lastg2 (this) |
Get the last g^2. | |
procedure | get_dg2 (this) |
Get the change in g^2. | |
procedure | eval_e (this, signals, a_model, gaussp, eq_steps, iou, eq_comm) |
Evaluate the e vector for a reconstruction step. | |
procedure | eval_f (this, derived_params, a_model) |
Evaluate the f vector for a reconstruction step. | |
procedure | eval_jacobians (this, signals, derived_params, locks, a_model, gaussp, params, eq_steps, iou, recon_comm, eq_comm) |
Evaluate the jacobian. | |
procedure | eval_step (this, signals, derived_params, locks, a_model, gaussp, params, eq_steps, iou, recon_comm, eq_comm) |
Evaluate the proper quasi-Newton step. | |
procedure | sl_step (this, k_use, step_size, delta_a) |
Take a straight line step. | |
procedure | seg_step (this, k_use, step_size, delta_a) |
Take a segmented path step. | |
procedure | lm_step (this, k_use, step_size, delta_a) |
Take a Levenberg-Marquardt step. | |
procedure | lm_rootfind (this, k_use, lambda_default, step_size, utdote) |
Root find for Levenberg-Marquardt step. | |
procedure | lm_function (this, f_sqrd, step_size_sqrd, lambda) |
Levenberg-Marquardt function. | |
procedure | step (this, signals, derived_params, locks, a_model, gaussp, params, eq_steps, iou, recon_comm, eq_comm) |
Step the reconstruction. | |
procedure | try_step (this, signals, derived_params, locks, a_model, gaussp, params, eq_steps, max_step, iou, recon_comm, eq_comm) |
Try step size. | |
procedure | eval_sem (this, params, signals, derived_params) |
Evaluate the parameter covariance and signal effectiveness. | |
procedure | invert_matrix (this, matrix, sub, name) |
Invert an M x N matrix. | |
procedure | write (this, iou) |
Write out a reconstruction. | |
procedure | write_step1 (this, step_type, iou) |
Write out a reconstruction step attempt. | |
procedure | write_step2 (this, iou) |
Write out a final reconstruction step. | |
GENERIC | write_step write_step1, write_step2 |
procedure | restart (this, result_ncid, current_step, signals, derived_params, a_model) |
Restart from result file. | |
procedure | sync_state (this, recon_comm) |
Syncronize state to children. | |
procedure | sync_svd (this, recon_comm) |
Syncronize svd. | |
procedure | sync_parent (this, params, num_signal, num_derived_params, stride, offset, recon_comm) |
Syncronize child state back to parent. | |
Public Attributes | |
integer | step_type = reconstruction_no_step_type |
Type descriptor of the boundry type for the lower(1) and upper(2) ranges. | |
real(rprec) | step_max |
Maximum number of reconstruction steps. | |
logical | use_central = .false. |
Controls if central differencing is used. | |
integer | current_step = 0 |
Counter to track the step number. | |
integer | last_para_signal = 0 |
Index of the last parallelable signal. | |
real(rprec), dimension(:,:), pointer | e => null() |
Array of all e vectors for each reconstruction step. | |
real(rprec), dimension(:,:), pointer | f => null() |
Array of all f vectors for each reconstruction step. The f vector contains the initial value of the derived parameter. | |
real(rprec), dimension(:,:), pointer | jacobian => null() |
The normalized jacobian for the current step. | |
real(rprec), dimension(:,:), pointer | derived_jacobian => null() |
The normalized derived jacobian for the current step. | |
real(rprec), dimension(:,:), pointer | hessian => null() |
The normalized hessian for the current step. | |
real(rprec), dimension(:), pointer | gradient => null() |
The normalized gradient for the current step. | |
real(rprec), dimension(:,:), pointer | delta_a => null() |
The normalized parameter step for each sigular value. | |
real(rprec), dimension(:), pointer | j_svd_w => null() |
The matrix of singular value decomposition singular values. | |
real(rprec), dimension(:,:), pointer | j_svd_u => null() |
The U matrix of singular value decomposition. | |
real(rprec), dimension(:,:), pointer | j_svd_vt => null() |
The tansposed V matrix of singular value decomposition. | |
real(rprec), dimension(:), pointer | delta_a_len => null() |
The normalized step length for each sigular value. | |
real(rprec), dimension(:), pointer | svd_w => null() |
The singular values used in inverting parameter covariance matrix. | |
real(rprec), dimension(:), pointer | exp_g2 => null() |
The expected g^2 for each reconstruction step. | |
real(rprec), dimension(:), pointer | step_size => null() |
The normalized step size for each reconstruction step. | |
integer, dimension(:), pointer | num_sv => null() |
The number of sigular values used for each reconstruction step. | |
real(rprec) | cut_svd |
Cutoff value for relative singular values. | |
real(rprec) | cut_eff |
Cutoff value for expected step efficiency. | |
real(rprec) | cut_marg_eff |
Cutoff value for expected marginal step efficiency. | |
real(rprec) | cut_delta_a |
Cutoff value for expected step size. | |
real(rprec) | cut_dg2 |
Cutoff value for expected change in g^2. | |
real(rprec) | cut_inv_svd |
Cutoff value for pseudo inverse singular values. | |
real(rprec), dimension(:,:), pointer | lm_ratio => null() |
Ratio of singular values to the Levenberg-Marquardt lambda. | |
real(rprec), dimension(:,:), pointer | last_values => null() |
Last signals. | |
Base class containing all the data needed to reconstruct a model.
procedure reconstruction::reconstruction_class::eval_e | ( | class (reconstruction_class), intent(inout) | this, |
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
class (model_class), pointer | a_model, | ||
type (gaussp_class_pointer), dimension(:), intent(inout) | gaussp, | ||
integer, intent(inout) | eq_steps, | ||
integer, intent(in) | iou, | ||
integer, intent(in) | eq_comm | ||
) |
Evaluate the e vector for a reconstruction step.
Evaluates the e vector after soving the equilibrium.
[in,out] | this | A reconstruction_class instance. |
[in,out] | signals | Array of signal objects. |
[in,out] | a_model | A model instance. |
[in,out] | gaussp | Array of guassian_process::gaussp_class objects. |
[in,out] | eq_steps | Number of steps for the equilibrium to take. |
[in] | iou | Input/output file to write log messages to. |
[in] | eq_comm | MPI communicator for the child equilibrium processes. |
procedure reconstruction::reconstruction_class::eval_f | ( | class (reconstruction_class), intent(inout) | this, |
type (param_pointer), dimension(:), intent(in) | derived_params, | ||
class (model_class), intent(inout) | a_model | ||
) |
Evaluate the f vector for a reconstruction step.
Evaluates the f vector after soving the equilibrium.
[in,out] | this | A reconstruction_class instance. |
[in] | derived_params | Array of param objects. |
[in,out] | a_model | A model instance. |
procedure reconstruction::reconstruction_class::eval_jacobians | ( | class (reconstruction_class), intent(inout) | this, |
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
type (param_pointer), dimension(:), intent(in) | derived_params, | ||
type (param_pointer), dimension(:), intent(in) | locks, | ||
class (model_class), pointer | a_model, | ||
type (gaussp_class_pointer), dimension(:), intent(inout) | gaussp, | ||
type (param_pointer), dimension(:), intent(inout) | params, | ||
integer, intent(inout) | eq_steps, | ||
integer, intent(in) | iou, | ||
integer, intent(in) | recon_comm, | ||
integer, intent(in) | eq_comm | ||
) |
Evaluate the jacobian.
Evaluates the jacobian. When using MPI, reach row of the jacobian is evaluated in parallel.
[in,out] | this | A reconstruction_class instance. |
[in,out] | signals | Array of signal objects. |
[in] | derived_params | Array of derived param objects. |
[in,out] | a_model | A model instance. |
[in,out] | params | Array of parameter objects. |
[in,out] | eq_steps | Number of steps for the equilibrium to take. |
[in] | iou | Input/output file to write log messages to. |
[in] | recon_comm | MPI communicator for the child jacobian processes. |
[in] | eq_comm | MPI communicator for the child equilibrium processes. |
procedure reconstruction::reconstruction_class::eval_sem | ( | class (reconstruction_class), intent(inout) | this, |
type (param_pointer), dimension(:), intent(inout) | params, | ||
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
type (param_pointer), dimension(:), intent(inout) | derived_params | ||
) |
Evaluate the parameter covariance and signal effectiveness.
Calculates the signal effectiveness matrix. This value is stored in parameter::param_class::sem. To find the signal effectiveness, an intermediate step computes the parameter covariance. Computes the derived parameter covariance matrix also.
[in,out] | this | A reconstruction_class instance. |
[in,out] | params | Array of parameter objects. |
[in,out] | signals | Array of signal objects. |
[in,out] | derived_params | Array of parameter objects. |
procedure reconstruction::reconstruction_class::eval_step | ( | class (reconstruction_class), intent(inout) | this, |
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
type (param_pointer), dimension(:), intent(in) | derived_params, | ||
type (param_pointer), dimension(:), intent(in) | locks, | ||
class (model_class), pointer | a_model, | ||
type (gaussp_class_pointer), dimension(:), intent(inout) | gaussp, | ||
type (param_pointer), dimension(:), intent(inout) | params, | ||
integer, intent(inout) | eq_steps, | ||
integer, intent(in) | iou, | ||
integer, intent(in) | recon_comm, | ||
integer, intent(in) | eq_comm | ||
) |
Evaluate the proper quasi-Newton step.
Evaluates the jacobian, performs the SVD and computes the quasi-Newton step.
[in,out] | this | A reconstruction_class instance. |
[in,out] | signals | Array of signal objects. |
[in] | derived_params | Array of derived param objects. |
[in,out] | a_model | A model instance. |
[in,out] | params | Array of parameter objects. |
[in,out] | eq_steps | Number of steps for the equilibrium to take. |
[in] | iou | Input/output file to write log messages to. |
[in] | recon_comm | MPI communicator for the child jacobian processes. |
[in] | eq_comm | MPI communicator for the child equilibrium processes. |
procedure reconstruction::reconstruction_class::get_dg2 | ( | class (reconstruction_class), intent(in) | this | ) |
Get the change in g^2.
Gets the change in g^2 from the last reconstruction step to the current step.
[in] | this | A reconstruction_class instance. |
procedure reconstruction::reconstruction_class::get_exp_dg2 | ( | class (reconstruction_class), intent(in) | this, |
real (rprec), dimension(:), intent(in) | delta_a | ||
) |
Get the expected dg^2 size.
Estimates the change in g^2. Equation 22 in Hanson et. al. doi:10.1088/0029-5515/49/7/075031 Note: there is a typo in equation 22. The - sign in the last term should be a plus sign.
[in] | this | A reconstruction_class instance. |
[in] | delta_a | The normalized parameter step. |
procedure reconstruction::reconstruction_class::get_exp_g2 | ( | class (reconstruction_class), intent(in) | this, |
real (rprec), dimension(:), intent(in) | delta_a | ||
) |
Get the expected g^2 size.
Estimates the g^2. Equation 22 in Hanson et. al. doi:10.1088/0029-5515/49/7/075031
[in] | this | A reconstruction_class instance. |
[in] | delta_a | The normalized parameter step. |
procedure reconstruction::reconstruction_class::get_g2 | ( | class (reconstruction_class), intent(in) | this | ) |
Get the current g^2.
Gets the g^2 value for the current reconstruction step.
[in] | this | A reconstruction_class instance. |
procedure reconstruction::reconstruction_class::get_k_use | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | num_sv | ||
) |
Get the number of singular values to use.
Number of sigular values to uses is determined by the various cutoffs.
[in,out] | this | A reconstruction_class instance. |
[in] | num_sv | Total number of sigular values. |
procedure reconstruction::reconstruction_class::get_lastg2 | ( | class (reconstruction_class), intent(in) | this | ) |
Get the last g^2.
Gets the g^2 value for the last reconstruction step.
[in] | this | A reconstruction_class instance. |
procedure reconstruction::reconstruction_class::invert_matrix | ( | class (reconstruction_class), intent(inout) | this, |
real (rprec), dimension(:,:), pointer | matrix, | ||
character (len=*), intent(in) | sub, | ||
character (len=*), intent(in) | name | ||
) |
Invert an M x N matrix.
Perform a pseudo-inversion using a singular value decomposition. This subroutine destroys the original matrix.
[in,out] | this | A reconstruction_class instance. |
[in,out] | matrix | The matrix to invert. |
[in] | sub | The subroutine of the matrix this was called from. |
[in] | name | The name of the matrix to be inverted. |
procedure reconstruction::reconstruction_class::lm_function | ( | class (reconstruction_class), intent(inout) | this, |
real (rprec), dimension(:), intent(in) | f_sqrd, | ||
real (rprec), intent(in) | step_size_sqrd, | ||
real (rprec), intent(in) | lambda | ||
) |
Levenberg-Marquardt function.
Function used by reconstruction_lm_rootfind for a Levenberg-Marquardt step.
[in,out] | this | A reconstruction_class instance. |
[in] | f_sqrd | f^2 |
[in] | step_size_sqrd | Square of the current step size. |
[in] | lambda | Lambda value. |
procedure reconstruction::reconstruction_class::lm_rootfind | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | k_use, | ||
real (rprec), intent(in) | lambda_default, | ||
real (rprec), intent(in) | step_size, | ||
real (rprec), dimension(:), intent(in) | utdote | ||
) |
Root find for Levenberg-Marquardt step.
Find the root of the reconstruction_lm_function for a Levenberg-Marquardt step.
[in,out] | this | A reconstruction_class instance. |
[in] | k_use | Number of singular values to use. |
[in] | lambda_default | The default lambda. |
[in] | step_size | Current step size limit. |
[in] | utdote | Vector of U^T * e |
procedure reconstruction::reconstruction_class::lm_step | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | k_use, | ||
real (rprec), intent(in) | step_size, | ||
real (rprec), dimension(:), intent(out) | delta_a | ||
) |
Take a Levenberg-Marquardt step.
Determines the normalized change in parameters for a Levenberg-Marquardt step.
[in,out] | this | A reconstruction_class instance. |
[in] | k_use | Number of singular values to use. |
[in] | step_size | Current step size limit. |
[out] | delta_a | Normalized change in parameter. |
|
final |
Deconstruct a reconstruction_class object.
Deallocates memory and uninitializes a reconstruction_class object.
[in,out] | this | A reconstruction_class instance. |
procedure reconstruction::reconstruction_class::restart | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | result_ncid, | ||
integer, intent(in) | current_step, | ||
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
type (param_pointer), dimension(:), intent(in) | derived_params, | ||
class (model_class), pointer | a_model | ||
) |
Restart from result file.
Restarts the reconstruction from the result file. Note, that there is an external step counter and an internal step counter. Internally we will restart current step at zero.
[in,out] | this | A reconstruction_class instance. |
[in] | result_ncid | NetCDF file id of the result file. |
[in] | current_step | Step index to read variables from. |
[in,out] | a_model | A model instance. |
procedure reconstruction::reconstruction_class::seg_step | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | k_use, | ||
real (rprec), intent(in) | step_size, | ||
real (rprec), dimension(:), intent(out) | delta_a | ||
) |
Take a segmented path step.
Determines the normalized change in parameters for a segmented path step.
[in,out] | this | A reconstruction_class instance. |
[in] | k_use | Number of singular values to use. |
[in] | step_size | Current step size limit. |
[out] | delta_a | Normalized change in parameter. |
procedure reconstruction::reconstruction_class::sl_step | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | k_use, | ||
real (rprec), intent(in) | step_size, | ||
real (rprec), dimension(:), intent(out) | delta_a | ||
) |
Take a straight line step.
Determines the normalized change in parameters for a straight line step.
[in,out] | this | A reconstruction_class instance. |
[in] | k_use | Number of singular values to use. |
[in] | step_size | Current step size limit. |
[out] | delta_a | Normalized change in parameter. |
procedure reconstruction::reconstruction_class::step | ( | class (reconstruction_class), intent(inout) | this, |
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
type (param_pointer), dimension(:), intent(in) | derived_params, | ||
type (param_pointer), dimension(:), intent(in) | locks, | ||
class (model_class), pointer | a_model, | ||
type (gaussp_class_pointer), dimension(:), intent(inout) | gaussp, | ||
type (param_pointer), dimension(:), intent(inout) | params, | ||
integer, intent(inout) | eq_steps, | ||
integer, intent(in) | iou, | ||
integer, intent(in) | recon_comm, | ||
integer, intent(in) | eq_comm | ||
) |
Step the reconstruction.
Advance the reconstruction by a single step. This is the main V3FIT algorithm.
[in,out] | this | A reconstruction_class instance. |
[in,out] | signals | Array of signal objects. |
[in,out] | locks | Array of lock objects. |
[in,out] | a_model | A model instance. |
[in,out] | gaussp | Array of guassian_process::gaussp_class objects. |
[in,out] | params | Array of v3fit_params::param_class objects. |
[in,out] | eq_steps | Number of steps the equilibrium took. |
[in] | iou | Input/output file to write log messages to. |
[in] | recon_comm | MPI communicator for the child jacobian processes. |
[in] | eq_comm | MPI communicator for the child equilibrium processes. |
procedure reconstruction::reconstruction_class::sync_parent | ( | class (reconstruction_class), intent(inout) | this, |
type (param_pointer), dimension(:), intent(inout) | params, | ||
integer, intent(in) | num_signal, | ||
integer, intent(in) | num_derived_params, | ||
integer, dimension(:), intent(in) | stride, | ||
integer, dimension(:), intent(in) | offset, | ||
integer, intent(in) | recon_comm | ||
) |
Syncronize child state back to parent.
Syncs data between the parent and child processes. If MPI support is not compiled in this subroutine reduces to a no op.
[in,out] | this | A reconstruction_class instance. |
[in] | num_params | Number of reconstruction parameters. |
[in] | num_signal | Number of signals. |
[in] | num_derived_params | Number of derived parameters. |
[in] | stride | Length work size for each process. |
[in] | offset | Start offset for each process. |
[in] | recon_comm | A MPI recon_comm handle. |
procedure reconstruction::reconstruction_class::sync_state | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | recon_comm | ||
) |
Syncronize state to children.
Syncs data between the parent and child processes. If MPI support is not compiled in this subroutine reduces to a no op.
[in,out] | this | A reconstruction_class instance. |
[in] | recon_comm | A MPI recon_comm handle. |
procedure reconstruction::reconstruction_class::sync_svd | ( | class (reconstruction_class), intent(inout) | this, |
integer, intent(in) | recon_comm | ||
) |
Syncronize svd.
Syncs the svd information to the child processes so these processes can sample parameter space for the step sizes.
[in,out] | this | A reconstruction_class instance. |
[in] | recon_comm | A MPI recon_comm handle. |
procedure reconstruction::reconstruction_class::try_step | ( | class (reconstruction_class), intent(inout) | this, |
type (signal_pointer), dimension(:), intent(inout) | signals, | ||
type (param_pointer), dimension(:), intent(in) | derived_params, | ||
type (param_pointer), dimension(:), intent(in) | locks, | ||
class (model_class), pointer | a_model, | ||
type (gaussp_class_pointer), dimension(:), intent(inout) | gaussp, | ||
type (param_pointer), dimension(:), intent(inout) | params, | ||
integer, intent(inout) | eq_steps, | ||
real (rprec), intent(inout) | max_step, | ||
integer, intent(in) | iou, | ||
integer, intent(in) | recon_comm, | ||
integer, intent(in) | eq_comm | ||
) |
Try step size.
Sample reconstruction by trying to take a specified step size.
[in,out] | this | A reconstruction_class instance. |
[in,out] | signals | Array of signal objects. |
[in,out] | locks | Array of lock objects. |
[in,out] | a_model | A model instance. |
[in,out] | params | Array of parameter objects. |
[in,out] | eq_steps | Number of steps for the equilibrium to take. |
[in] | max_step | Upper limit of the step size to attemp. |
[in] | iou | Input/output file to write log messages to. |
[in] | recon_comm | MPI communicator for the child jacobian processes. |
[in] | eq_comm | MPI communicator for the child equilibrium processes. |
[out] | param_value | Parameter corresponding to the smallest g^2 value. |
procedure reconstruction::reconstruction_class::write | ( | class (reconstruction_class), intent(in) | this, |
integer, intent(in) | iou | ||
) |
Write out a reconstruction.
Write out the final information about a reconstruction. This writes out a table of the g^2, expected g^2, number of sigular values used and the normalized step size used.
[in] | this | A reconstruction_class instance. |
[in] | iou | Input/output file to write final result to. |
procedure reconstruction::reconstruction_class::write_step1 | ( | class (reconstruction_class), intent(in) | this, |
character (len=*), intent(in) | step_type, | ||
integer, intent(in) | iou | ||
) |
Write out a reconstruction step attempt.
Write out information about a reconstruction step. This writes the number of singular values use, the step type, the normalized step size and the expected g^2.
[in] | this | A reconstruction_class instance. |
[in] | step_type | The step type used. |
[in] | iou | Input/output file to write final result to. |
procedure reconstruction::reconstruction_class::write_step2 | ( | class (reconstruction_class), intent(in) | this, |
integer, intent(in) | iou | ||
) |
Write out a final reconstruction step.
Write out information about a reconstruction step. This writes the number of singular values use, the step type, the normalized step size and the expected g^2.
[in] | this | A reconstruction_class instance. |
[in] | iou | Input/output file to write final result to. |
real (rprec) reconstruction::reconstruction_class::cut_delta_a |
Cutoff value for expected step size.
real (rprec) reconstruction::reconstruction_class::cut_dg2 |
Cutoff value for expected change in g^2.
real (rprec) reconstruction::reconstruction_class::cut_eff |
Cutoff value for expected step efficiency.
real (rprec) reconstruction::reconstruction_class::cut_inv_svd |
Cutoff value for pseudo inverse singular values.
real (rprec) reconstruction::reconstruction_class::cut_marg_eff |
Cutoff value for expected marginal step efficiency.
real (rprec) reconstruction::reconstruction_class::cut_svd |
Cutoff value for relative singular values.
real (rprec), dimension(:,:), pointer reconstruction::reconstruction_class::e => null() |
Array of all e vectors for each reconstruction step.
real (rprec) reconstruction::reconstruction_class::step_max |
Maximum number of reconstruction steps.
integer reconstruction::reconstruction_class::step_type = reconstruction_no_step_type |
Type descriptor of the boundry type for the lower(1) and upper(2) ranges.