halo_module Module

Main module with all the stuff to solve the Halo problem.

main program


Uses

  • module~~halo_module~~UsesGraph module~halo_module halo_module bspline_module bspline_module module~halo_module->bspline_module ddeabm_module ddeabm_module module~halo_module->ddeabm_module fortran_astrodynamics_toolkit fortran_astrodynamics_toolkit module~halo_module->fortran_astrodynamics_toolkit json_module json_module module~halo_module->json_module module~config_file_module config_file_module module~halo_module->module~config_file_module module~halo_utilities_module halo_utilities_module module~halo_module->module~halo_utilities_module module~parameters_module parameters_module module~halo_module->module~parameters_module module~splined_ephemeris_module splined_ephemeris_module module~halo_module->module~splined_ephemeris_module nlesolver_module nlesolver_module module~halo_module->nlesolver_module numerical_differentiation_module numerical_differentiation_module module~halo_module->numerical_differentiation_module pyplot_module pyplot_module module~halo_module->pyplot_module module~config_file_module->json_module module~config_file_module->module~parameters_module module~halo_utilities_module->fortran_astrodynamics_toolkit module~halo_utilities_module->module~parameters_module module~halo_utilities_module->module~splined_ephemeris_module iso_fortran_env iso_fortran_env module~halo_utilities_module->iso_fortran_env module~halo_kinds_module halo_kinds_module module~halo_utilities_module->module~halo_kinds_module module~parameters_module->fortran_astrodynamics_toolkit module~parameters_module->module~halo_kinds_module module~splined_ephemeris_module->bspline_module module~splined_ephemeris_module->fortran_astrodynamics_toolkit module~splined_ephemeris_module->iso_fortran_env module~splined_ephemeris_module->module~halo_kinds_module module~halo_kinds_module->iso_fortran_env

Derived Types

type, public ::  segment_data

data for a segment

Components

Type Visibility Attributes Name Initial
real(kind=wp), public :: et_ref = zero

ephemeris time reference epoch (all times relative to this)

real(kind=wp), public :: t0 = zero

initial time [days]

real(kind=wp), public :: t0_scale = one

opt var scale for t0

real(kind=wp), public :: t0_dpert = 1.0e-5_wp

opt var dpert for t0 [scaled]

real(kind=wp), public, dimension(6) :: x0_rotating = zero

initial state [km, km/s], rotating frame

real(kind=wp), public, dimension(6) :: x0_rotating_scale = one

opt var scale for x0_rotating

real(kind=wp), public, dimension(6) :: x0_rotating_dpert = [1.0e-2_wp, 1.0e-2_wp, 1.0e-2_wp, 1.0e-5_wp, 1.0e-5_wp, 1.0e-5_wp]

opt var dperts for x0_rotating [scaled]

real(kind=wp), public, dimension(6) :: x0 = zero

initial state [km, km/s], inertial frame

real(kind=wp), public :: tf = zero

final time [days]

real(kind=wp), public, dimension(6) :: xf = zero

final state [km, km/s], inertial frame

real(kind=wp), public, dimension(6) :: xf_rotating = zero

final state [km, km/s], rotating frame

real(kind=wp), public, dimension(6) :: xf_rotating_scale = one

func scale for xf_rotating

type, public ::  trajectory

a segment trajectory (for plotting or output)

Components

Type Visibility Attributes Name Initial
real(kind=wp), public, dimension(:), allocatable :: et
real(kind=wp), public, dimension(:), allocatable :: x
real(kind=wp), public, dimension(:), allocatable :: y
real(kind=wp), public, dimension(:), allocatable :: z
real(kind=wp), public, dimension(:), allocatable :: vx
real(kind=wp), public, dimension(:), allocatable :: vy
real(kind=wp), public, dimension(:), allocatable :: vz

Type-Bound Procedures

procedure, public :: destroy => destroy_trajectory
procedure, public :: add => add_point_to_trajectory

type, public, extends(ddeabm_class) ::  segment

a ballistic segment in the mission

Read more…

Components

Type Visibility Attributes Name Initial
character(len=:), public, allocatable :: name

the segment name

type(segment_data), public :: data
type(segment_data), public :: cached_data

used when computing gradients

type(geopotential_model_pines), public, pointer :: grav => null()

central body geopotential model [global]

class(jpl_ephemeris), public, pointer :: eph => null()

the ephemeris [global]

logical, public :: pointmass_central_body = .false.
type(trajectory), public :: traj_inertial

in the inertial frame (j2000-moon)

type(trajectory), public :: traj_rotating

in the rotating frame (moon-earth, moon-centered)

type(trajectory), public :: traj_se_rotating

in the rotating frame (sun-earth, earth-centered)

Type-Bound Procedures

procedure, public :: set_input => set_segment_inputs
procedure, public :: get_inputs => get_segment_inputs
procedure, public :: get_outputs => get_segment_outputs
procedure, public :: set_outputs => set_segment_outputs
procedure, public :: propagate => propagate_segment
procedure, public :: cache
procedure, public :: uncache
procedure, public :: put_data_in_segment

type, public, extends(numdiff_type) ::  mission_type

the mission [this is a numdiff_type for organizational purposes.... rethink this maybe ...

Components

Type Visibility Attributes Name Initial
character(len=:), public, allocatable :: initial_guess_from_file

to read the initial guess from a JSOn file

logical, public :: use_splined_ephemeris = .false.

if true, the ephemeris is splined

real(kind=wp), public :: dt_spline_sec = 3600.0_wp

time step in seconds for spline step [1 hr]

logical, public :: pointmass_central_body = .false.

if true, the central body is a pointmass (moon). otherwise, the grav model is used.

type(geopotential_model_pines), public, pointer :: grav => null()

central body geopotential model [global]

class(jpl_ephemeris), public, pointer :: eph => null()

the ephemeris [global]

type(segment), public, dimension(:), allocatable :: segs

the list of segments

real(kind=wp), public, dimension(:), allocatable :: xscale

opt var scale factors

real(kind=wp), public, dimension(:), allocatable :: dpert_

opt var dpert factors [rename because already in base class]

real(kind=wp), public, dimension(:), allocatable :: fscale

function scale factors

character(len=20), public, dimension(:), allocatable :: xname

opt var names

character(len=:), public, allocatable :: N_or_S

'N' or 'S'

character(len=:), public, allocatable :: L1_or_L2

'L1' or 'L2'

logical, public :: solve = .true.

run the nle solver to solve the problem

logical, public :: generate_plots = .true.

to generate the python plots

logical, public :: generate_trajectory_files = .true.

to export the txt trajectory files.

logical, public :: generate_json_trajectory_file = .false.

to export the JSON trajectory file for the solution (contains inertial and rotating data). This one is used by the PyVista python script for plotting

logical, public :: generate_guess_and_solution_files = .true.

to export the json guess and solution files.

logical, public :: generate_kernel = .false.

to generate a spice kernel (bsp) of the solution [this requires the external mkspk SPICE tool]

logical, public :: generate_defect_file = .false.

generate the file that shows the pos/vel constraint defects for the solution

logical, public :: generate_eclipse_files = .false.

generate the eclipse data file for the solution

logical, public :: run_pyvista_script = .false.

run the pyvista script to generate the interacdtive 3d plot

real(kind=wp), public :: r_eclipse_bubble = 0.0_wp

radius of the "eclipse bubble" [km]

real(kind=wp), public :: eclipse_dt_step = 3600.0_wp

dense time step output for eclipse calculations in sec [1 hour]

integer, public :: eclipse_filetype = 1

type of eclipse file: 1=csv or 2=json

integer, public :: epoch_mode = 1

1 - calendar date was specified, 2 - ephemeris time was specified

integer, public :: year = 0

epoch of first point (first periapsis crossing)

integer, public :: month = 0

epoch of first point (first periapsis crossing)

integer, public :: day = 0

epoch of first point (first periapsis crossing)

integer, public :: hour = 0

epoch of first point (first periapsis crossing)

integer, public :: minute = 0

epoch of first point (first periapsis crossing)

real(kind=wp), public :: sec = zero

epoch of first point (first periapsis crossing)

real(kind=wp), public :: et_ref = zero

ephemeris time reference epoch [sec] (all times relative to this). this is computed from the year,month,day,hour,minute,sec values if epoch_mode=1

integer, public :: n_revs = 10

Number of revs in the Halo.

real(kind=wp), public :: rtol = 1.0e-12_wp

integrator rtol

real(kind=wp), public :: atol = 1.0e-12_wp

integrator atol

real(kind=wp), public :: nlesolver_tol = 1.0e-6_wp

nlesolver tol

real(kind=wp), public, dimension(:), allocatable :: initial_r

if fix_initial_r is True, this can be used to specify the initial r (moon-centered earth-moon rotating frame)

character(len=:), public, allocatable :: ephemeris_file

the JPL ephemeris file to load

character(len=:), public, allocatable :: gravfile

spherical harmonic gravity coeff file (Moon) example: 'data/grav/gggrx_0020pm_sha.tab'

character(len=:), public, allocatable :: patch_point_file

Halo CR3BP patch point solution file example: 'data/L2_halos.json'

logical, public :: patch_point_file_is_periapsis = .false.

if the state in the patch point file are periapsis states (false means they are apoapsis states. the ones in L2_halos.json are apoapsis states)

type(patch_point), public :: periapsis

Patchpoint Periapse State

type(patch_point), public :: quarter

Patchpoint 1/4 Rev State

type(patch_point), public :: apoapsis

Patchpoint 1/2 Rev State

real(kind=wp), public :: period = 0.0_wp

Halo period [days]

real(kind=wp), public :: period8 = 0.0_wp

1/8 of Halo period [days]

logical, public :: fix_initial_time = .false.

to fix the initial epoch in the mission

logical, public :: fix_initial_r = .false.

fix the initial periapsis position vector (x,y,z)

integer, public :: fix_ry_at_end_of_rev = 0

fix ry at the end of the specified rev (at periapsis) not used if <= 0.

logical, public :: fix_final_ry_and_vx = .false.

fix ry and vx at the end of the last rev (at periapsis)

integer, public :: solver_mode = 1

use sparse or dense solver:

Read more…

Type-Bound Procedures

procedure, public :: init => initialize_the_mission
procedure, public :: plot => plot_trajectory
procedure, public :: export_trajectory_json_file
procedure, public :: write_optvars_to_file
procedure, public :: constraint_violations
procedure, public :: print_constraint_defects
procedure, public :: generate_eclipse_data
procedure, public :: get_problem_arrays
procedure, public :: put_x_in_segments
procedure, public :: get_x_from_json_file
procedure, public :: get_scales_from_segs
procedure, public :: segs_to_propagate
procedure, public :: get_sparsity_pattern
procedure, public :: define_problem_size
procedure, public :: get_case_name
procedure, public :: generate_patch_points
procedure, public :: update_epoch

type, public, extends(nlesolver_type) ::  my_solver_type

the solver class, which contains an instance of the mission

Components

Type Visibility Attributes Name Initial
type(mission_type), public :: mission

the Halo mission

Type-Bound Procedures

procedure, public :: init => initialize_the_solver
procedure, public :: read_config_file

Functions

public function get_case_name(me) result(case_name)

Returns a string that can be used for file names, etc. (read the values of the global variables).

Read more…

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(in) :: me

Return Value character(len=:), allocatable


Subroutines

public subroutine halo_solver_main(config_file_name, debug)

Main program to solve the Halo targeting problem.

Read more…

Arguments

Type IntentOptional Attributes Name
character(len=*), intent(in) :: config_file_name

the config file to read

logical, intent(in) :: debug

for debugging prints

public pure subroutine destroy_trajectory(me)

Destroy a trajectory (deallocate the data arrays).

Arguments

Type IntentOptional Attributes Name
class(trajectory), intent(inout) :: me

public pure subroutine define_problem_size(me, n, m, n_segs, n_nonzero, full_problem)

Returns the size variables for the "forward-backward" formulation of the problem.

Read more…

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(in) :: me
integer, intent(out), optional :: n

number of opt vars for the solver

integer, intent(out), optional :: m

number of equality constraints for the solver

integer, intent(out), optional :: n_segs

number of segments

integer, intent(out), optional :: n_nonzero

number of nonzero elements in the Jacobian

logical, intent(in), optional :: full_problem

if true (default) always return the results for the full problem (i.e., without fixing any of the variables)

public subroutine put_data_in_segment(me, d)

Set all the data in a segment.

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me
type(segment_data), intent(in) :: d

public subroutine cache(me)

Cache all the data in a segment.

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me

public subroutine uncache(me)

Restore all the segment data form the cache.

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me

public subroutine initialize_the_solver(me, config_file_name, x)

Initialize the mission.

Read more…

Arguments

Type IntentOptional Attributes Name
class(my_solver_type), intent(inout) :: me
character(len=*), intent(in) :: config_file_name

the config file to read

real(kind=wp), intent(out), optional, dimension(:), allocatable :: x

initial guess

public subroutine qrm_solver(me, n_cols, n_rows, n_nonzero, irow, icol, a, b, x, istat)

Custom solver that uses QR_MUMPS

Read more…

Arguments

Type IntentOptional Attributes Name
class(nlesolver_type), intent(inout) :: me
integer, intent(in) :: n_cols

n: number of columns in A.

integer, intent(in) :: n_rows

m: number of rows in A.

integer, intent(in) :: n_nonzero

number of nonzero elements of A.

integer, intent(in), dimension(n_nonzero) :: irow

sparsity pattern (size is n_nonzero)

integer, intent(in), dimension(n_nonzero) :: icol

sparsity pattern (size is n_nonzero)

real(kind=wp), intent(in), dimension(n_nonzero) :: a

matrix elements (size is n_nonzero)

real(kind=wp), intent(in), dimension(n_rows) :: b

right hand side (size is m)

real(kind=wp), intent(out), dimension(n_cols) :: x

solution (size is n)

integer, intent(out) :: istat

status code (=0 for success)

public subroutine get_problem_arrays(me, x, f)

Get the arrays for the problem

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(in) :: me
real(kind=wp), intent(out), optional, dimension(:) :: x

opt var vector [scaled]

real(kind=wp), intent(out), optional, dimension(:) :: f

constraint violations

public subroutine print_constraint_defects(me, filename)

Print the r and v constraint defect norms for each segment constraint.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(in) :: me
character(len=*), intent(in) :: filename

csv file to write to

public subroutine get_x_from_json_file(me, x)

Read the JSON solution file and put the x vector in the mission. Can be used to restart a solution from a previous run (e.g., with different settings, as long as the fundamental problem isn't changed).

Read more…

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
real(kind=wp), intent(out), dimension(:), allocatable :: x

scaled x vector

public subroutine put_x_in_segments(me, x_scaled)

Take the x vector from the solver, and use it to populate all the segments

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x_scaled

opt var vector (scaled)

public subroutine get_sparsity_pattern(me, irow, icol, linear_irow, linear_icol, linear_vals, maxgrp, ngrp)

Returns the sparsity pattern for the "forward-backward" Halo problem.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
integer, intent(out), dimension(:), allocatable :: irow

sparsity pattern nonzero elements row indices

integer, intent(out), dimension(:), allocatable :: icol

sparsity pattern nonzero elements column indices

integer, intent(out), optional, dimension(:), allocatable :: linear_irow

linear sparsity pattern nonzero elements row indices

integer, intent(out), optional, dimension(:), allocatable :: linear_icol

linear sparsity pattern nonzero elements column indices

real(kind=wp), intent(out), optional, dimension(:), allocatable :: linear_vals

linear sparsity values (constant elements of the Jacobian)

integer, intent(out), optional :: maxgrp

DSM sparsity partition

integer, intent(out), optional, dimension(:), allocatable :: ngrp

DSM sparsity partition

public subroutine initialize_the_mission(me, x)

Initialize the mission.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
real(kind=wp), intent(out), optional, dimension(:), allocatable :: x

initial guess

public subroutine get_scales_from_segs(me)

Populate the xscale and fscale problem arrays from the segment data.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me

public subroutine ballistic_derivs(me, t, x, xdot)

Equations of motion for a ballistic segment.

Read more…

Arguments

Type IntentOptional Attributes Name
class(ddeabm_class), intent(inout) :: me
real(kind=wp), intent(in) :: t

time [sec from epoch]

real(kind=wp), intent(in), dimension(:) :: x

state [r,v] in inertial frame (moon-centered)

real(kind=wp), intent(out), dimension(:) :: xdot

derivative of state ()

public subroutine set_segment_inputs(me, t0, tf, x0_rotating)

Sets all the info in a segment for it to be propagated.

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me
real(kind=wp), intent(in) :: t0

[days]

real(kind=wp), intent(in) :: tf

[days]

real(kind=wp), intent(in), dimension(6) :: x0_rotating

state in rotating frame

public subroutine get_segment_inputs(me, t0, x0_rotating)

Gets the initial states of a segment

Arguments

Type IntentOptional Attributes Name
class(segment), intent(in) :: me
real(kind=wp), intent(out), optional :: t0
real(kind=wp), intent(out), optional, dimension(6) :: x0_rotating

rotating frame

public subroutine get_segment_outputs(me, xf, xf_rotating, x0_rotating)

After propagating a segment, this gets the outputs.

Arguments

Type IntentOptional Attributes Name
class(segment), intent(in) :: me
real(kind=wp), intent(out), optional, dimension(6) :: xf

inertial frame

real(kind=wp), intent(out), optional, dimension(6) :: xf_rotating

rotating frame

real(kind=wp), intent(out), optional, dimension(6) :: x0_rotating

public subroutine set_segment_outputs(me, xf, xf_rotating)

Set the outputs of a segment, assuming it has been propagated elsewhere

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me
real(kind=wp), intent(in), dimension(6) :: xf

inertial frame

real(kind=wp), intent(in), dimension(6) :: xf_rotating

inertial frame

public subroutine propagate_segment(me, mode, tstep)

Propagate a segment (assumes the inputs have already been populated)

Arguments

Type IntentOptional Attributes Name
class(segment), intent(inout) :: me
integer, intent(in), optional :: mode

1 - don't report steps, 2 - report steps (for plotting)

real(kind=wp), intent(in), optional :: tstep

fixed time step for mode=2

public subroutine segs_to_propagate(me, funcs_to_compute, isegs_to_propagate)

Given the functions to be computed, returns the segments that need to be propagated.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
integer, intent(in), optional, dimension(:) :: funcs_to_compute

the indices of f to compute

integer, intent(out), dimension(:), allocatable :: isegs_to_propagate

segment numbers

public subroutine constraint_violations(me, x, f, funcs_to_compute)

Computes the constraint violation vector for the mission.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x

opt var vector for the mission [n]

real(kind=wp), intent(out), dimension(:) :: f

constraint violation vector for the mission [m]

integer, intent(in), optional, dimension(:) :: funcs_to_compute

the indices of f to compute

public subroutine halo_func(me, x, f)

Compute the solver function (all the constraint violations)

Arguments

Type IntentOptional Attributes Name
class(nlesolver_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(out), dimension(:) :: f

public subroutine halo_grad(me, x, g)

Compute the gradient of the solver function (Jacobian matrix). Dense version.

Arguments

Type IntentOptional Attributes Name
class(nlesolver_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(out), dimension(:,:) :: g

public subroutine halo_grad_sparse(me, x, g)

Compute the gradient of the solver function (Jacobian matrix). Sparse version.

Arguments

Type IntentOptional Attributes Name
class(nlesolver_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(out), dimension(:) :: g

public subroutine halo_grad_info(me, column, i, x)

Info function for the gradients.

Read more…

Arguments

Type IntentOptional Attributes Name
class(numdiff_type), intent(inout) :: me
integer, intent(in), dimension(:) :: column

the columns being computed.

integer, intent(in) :: i

perturbing these columns for the ith time (1,2,...)

real(kind=wp), intent(in), dimension(:) :: x

the nominal variable vector

public subroutine my_func(me, x, f, funcs_to_compute)

Wrapper for the constraint violation function (use by NumDiff)

Arguments

Type IntentOptional Attributes Name
class(numdiff_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(out), dimension(:) :: f
integer, intent(in), dimension(:) :: funcs_to_compute

the function indices to compute in the full f vector

public subroutine halo_export(me, x, f, iter)

export an iteration from the solver

Arguments

Type IntentOptional Attributes Name
class(nlesolver_type), intent(inout) :: me
real(kind=wp), intent(in), dimension(:) :: x
real(kind=wp), intent(in), dimension(:) :: f
integer, intent(in) :: iter

iteration number

public subroutine plot_trajectory(me, filename, export_trajectory, draw_trajectory, only_first_rev)

Plot the trajectory using matplotlib and/or generate a text file (for MKSPK)

Read more…

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
character(len=*), intent(in) :: filename

plot file name [without extension]

logical, intent(in), optional :: export_trajectory

if true, a text trajectory file is also produced (that can be read by MKSPK) [default is False]

logical, intent(in), optional :: draw_trajectory

if true [Default], a matplotlib plot is generated

logical, intent(in), optional :: only_first_rev

to only do the first rev [Default is False]

public subroutine export_trajectory_json_file(me, filename, only_first_rev)

Export the trajectory JSON file.

Read more…

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
character(len=*), intent(in) :: filename

plot file name [without extension]

logical, intent(in), optional :: only_first_rev

to only do the first rev [Default is False]

public subroutine generate_eclipse_data(me, fileprefix, filetype)

Propagate all the segments with a dense time step and export the Eclipsing data w.r.t. Earth. In this file, any point <0 is in eclipse.

Read more…

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
character(len=*), intent(in) :: fileprefix

file prefix for the csv file (case name will be added)

integer, intent(in), optional :: filetype

type of output file: 1=csv [default] or 2=json

public subroutine add_point_to_trajectory(me, et, x)

Add a point to a trajectory variable.

Arguments

Type IntentOptional Attributes Name
class(trajectory), intent(inout) :: me

the trajectory

real(kind=wp), intent(in) :: et
real(kind=wp), intent(in), dimension(6) :: x

public subroutine trajectory_export_func(me, t, x)

Export a trajectory point during propagation.

Read more…

Arguments

Type IntentOptional Attributes Name
class(ddeabm_class), intent(inout) :: me
real(kind=wp), intent(in) :: t

time from integrator (sec from sim epoch)

real(kind=wp), intent(in), dimension(:) :: x

state [moon-centered inertial frame]

public subroutine read_epoch(f, epoch_mode, year, month, day, hour, minute, sec, et)

Read the epoch info from a config file.

Arguments

Type IntentOptional Attributes Name
type(config_file), intent(inout) :: f
integer, intent(out) :: epoch_mode

1 : calendar date was specified, 2: ephemeris time was specified

integer, intent(out) :: year

only output if found (epoch_mode=1)

integer, intent(out) :: month

only output if found (epoch_mode=1)

integer, intent(out) :: day

only output if found (epoch_mode=1)

integer, intent(out) :: hour

only output if found (epoch_mode=1)

integer, intent(out) :: minute

only output if found (epoch_mode=1)

real(kind=wp), intent(out) :: sec

only output if found (epoch_mode=1)

real(kind=wp), intent(out) :: et

the ephemeris time [always output]

public subroutine read_config_file(me, filename)

Read the config file that defines the problem to be solved and set all the global variables.

Read more…

Arguments

Type IntentOptional Attributes Name
class(my_solver_type), intent(inout) :: me
character(len=*), intent(in) :: filename

the JSON config file to read

public subroutine update_epoch(me)

Update the variables for the reference epoch in the mission class. Either computes the et from the calendar state, or vice versa.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me

public subroutine generate_patch_points(me, lstar, tstar, period, x0, z0, ydot0, pp)

Generates the patch points (unnormalized, moon-centered) from the CR3BP normalized guess.

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
real(kind=wp), intent(in) :: lstar
real(kind=wp), intent(in) :: tstar
real(kind=wp), intent(in) :: period
real(kind=wp), intent(in) :: x0
real(kind=wp), intent(in) :: z0
real(kind=wp), intent(in) :: ydot0
type(patch_point), intent(out), dimension(3) :: pp

periapsis, quarter, apoapsis patch points

public subroutine write_optvars_to_file(me, filename, x)

Write the current x variables to a file (compatible with the Copernicus optvar file).

Arguments

Type IntentOptional Attributes Name
class(mission_type), intent(inout) :: me
character(len=*), intent(in) :: filename

file name [without extension]

real(kind=wp), intent(in), dimension(:) :: x

solver opt vars vector