Module MAPLEAF.IO

Input/Output functionality:

  • Reading/Writing Simulation Definition Files
  • Plotting/Logging results

MAPLEAF.IO Relies on MAPLEAF.Motion to implement a few convenience / parsing functions.

Expand source code
'''
Input/Output functionality:

* Reading/Writing Simulation Definition Files
* Plotting/Logging results

MAPLEAF.IO Relies on MAPLEAF.Motion to implement a few convenience / parsing functions.

.. image:: https://upload.wikimedia.org/wikipedia/commons/thumb/f/f8/Laptop-hard-drive-exposed.jpg/1280px-Laptop-hard-drive-exposed.jpg
'''
# Make the classes in all submodules importable directly from MAPLEAF.Rocket
from .rocketFlight import *
from .simDefinition import *
from .subDictReader import *
from .CythonLog import *

subModules = [ rocketFlight, simDefinition, subDictReader, CythonLog ]

__all__ = []

for subModule in subModules:
    __all__ += subModule.__all__

Sub-modules

MAPLEAF.IO.CythonLog
MAPLEAF.IO.HIL

Hardware in the loop simulation functionality. Orientation/Position/Velocity info is passed to a hardware control system, which passes control inputs …

MAPLEAF.IO.Logging

Classes and functions for creating simulation logs for regular simulations (Logger) and Monte Carlo simulations (MonteCarloLogger)

MAPLEAF.IO.Plotting

Functions to create plots and animations of simulations results.

MAPLEAF.IO.gridConvergenceFunctions

Functions to analyze grid convergence of simulations

MAPLEAF.IO.rocketFlight

Temporarily holds simulation results. Used for creating flight animations.

MAPLEAF.IO.simDefinition

Contains a class to read, write and modify simulation definition (.mapleaf) files, the master dictionary of default values for simulation …

MAPLEAF.IO.subDictReader

Wrapper class to read from a specific sub-dictionary in a SimDefinition.

Functions

def getAbsoluteFilePath(relativePath: str, alternateRelativeLocation: str = '', silent=False) ‑> str

Takes a path defined relative to the MAPLEAF repository and tries to return an absolute path for the current installation. alternateRelativeLocation (str) location of an alternate file/folder the path could be relative to Returns original relativePath if an absolute path is not found

Expand source code
def getAbsoluteFilePath(relativePath: str, alternateRelativeLocation: str = "", silent=False) -> str:
    ''' 
        Takes a path defined relative to the MAPLEAF repository and tries to return an absolute path for the current installation.
        alternateRelativeLocation (str) location of an alternate file/folder the path could be relative to
        Returns original relativePath if an absolute path is not found
    '''
    # Check if path is relative to MAPLEAF installation location
    # This file is at MAPLEAF/IO/SimDefinition, so MAPLEAF's install directory is three levels up
    pathToMAPLEAFInstallation = Path(__file__).parent.parent.parent

    relativePath = Path(relativePath)
    absolutePath = pathToMAPLEAFInstallation / relativePath

    if absolutePath.exists():
        return str(absolutePath)
    else:
        if alternateRelativeLocation != "":
            # Try alternate location
            alternateLocation = Path(alternateRelativeLocation)
            
            if alternateLocation.is_file():
                # If the alternate location provided is a file path, check if the file is in the parent directory
                absolutePath = alternateLocation.parent / relativePath
            else:
                absolutePath = alternateLocation / relativePath

            if absolutePath.exists():
                return str(absolutePath)
                
    if not silent:
        if ".pdf" not in str(relativePath): # Assume .pdf files are outputs, may not be created yet - so we wouldn't expect to find them immediately!
            print("WARNING: Unable to compute absolute path replacement for: {}, try providing an absolute path".format(relativePath))
        
    return str(relativePath)

Classes

class Log (...)

Class manages logs for any number of states that need to be logged. Each state/parameter that needs to be logged is given its own column in the log (represented as a python list) Each row holds values from a single time. When a new row is added to the log, any empty values in the previous row will be filled with the fill value.

Internally, the log is represented as a dict containing a list for each log column

columnNames should be a list of strings To get references to the log columns, initialize the log as empty and then use the addColumn(s) functions. Those return the list(s) they add to the log for direct access

Subclasses

Instance variables

var fillValues

Return an attribute of instance, which is of type owner.

var logColumns

Return an attribute of instance, which is of type owner.

Methods

def addColumn(...)

Returns reference to the log column (list)

def addColumns(...)

Returns list of references to the log columns (list[list])

def deleteLastRow(...)
def expandIterableColumns(...)

Look for columns containing vector values, convert them into a triplet of values

def getValue(...)
def logValue(...)

Throws KeyError if colName is not a valid column name

def newLogRow(...)

Start a new row in the log. Typically called after a time step or force evaluation is completed

def writeToCSV(...)
class RocketFlight

Holds simulation results

These arrays are filled in during a flight simulation. If the rocket is not controlled, self.actuatorDefls will be set to None

Expand source code
class RocketFlight():
    ''' Holds simulation results '''

    def __init__(self):
        '''  These arrays are filled in during a flight simulation. If the rocket is not controlled, self.actuatorDefls will be set to None '''
        self.times = []
        self.rigidBodyStates = []
        self.actuatorDefls = []
        self.actuatorTargetDefls = []

        # To be filled during flight
        self.engineOffTime = None
        self.mainChuteDeployTime = None
        self.targetLocation = None

    def getFlightTime(self):
        return self.times[-1]

    def getApogee(self):
        maxAltitude = -10000000
        for state in self.rigidBodyStates:
            if state.position.Z > maxAltitude:
                maxAltitude = state.position.Z
        return maxAltitude

    def getMaxSpeed(self):
        maxVel = 0
        for state in self.rigidBodyStates:
            if state.velocity.length() > maxVel:
                maxVel = state.velocity.length()
        return maxVel

    def getMaxHorizontalVel(self):
        maxHorizontalVel = 0
        for state in self.rigidBodyStates:
            horizVel = sqrt(state.velocity.X**2 + state.velocity.Y**2)
            if horizVel > maxHorizontalVel:
                maxHorizontalVel = horizVel
        return maxHorizontalVel

    def getLandingLocation(self):
        '''
            Extrapolates to z = 0
            If the simulation stopped before z = 0 is reached, the last position is returned
        '''

        #TODO: Handle non-flat surface to land on
        secondLastPosition = self.rigidBodyStates[-2].position
        lastPosition = self.rigidBodyStates[-1].position

        dz = lastPosition.Z - secondLastPosition.Z

        if lastPosition.Z > dz*2:
            dxdz = (lastPosition.X - secondLastPosition.X) / dz
            dydz = (lastPosition.Y - secondLastPosition.Y) / dz

            dzLand = -lastPosition.Z 

            xLand = lastPosition.X + dxdz*dzLand
            yLand = lastPosition.Y + dydz*dzLand
            zLand = 0.0
            
            return Vector(xLand, yLand, zLand)
        else:
            return lastPosition

Methods

def getApogee(self)
Expand source code
def getApogee(self):
    maxAltitude = -10000000
    for state in self.rigidBodyStates:
        if state.position.Z > maxAltitude:
            maxAltitude = state.position.Z
    return maxAltitude
def getFlightTime(self)
Expand source code
def getFlightTime(self):
    return self.times[-1]
def getLandingLocation(self)

Extrapolates to z = 0 If the simulation stopped before z = 0 is reached, the last position is returned

Expand source code
def getLandingLocation(self):
    '''
        Extrapolates to z = 0
        If the simulation stopped before z = 0 is reached, the last position is returned
    '''

    #TODO: Handle non-flat surface to land on
    secondLastPosition = self.rigidBodyStates[-2].position
    lastPosition = self.rigidBodyStates[-1].position

    dz = lastPosition.Z - secondLastPosition.Z

    if lastPosition.Z > dz*2:
        dxdz = (lastPosition.X - secondLastPosition.X) / dz
        dydz = (lastPosition.Y - secondLastPosition.Y) / dz

        dzLand = -lastPosition.Z 

        xLand = lastPosition.X + dxdz*dzLand
        yLand = lastPosition.Y + dydz*dzLand
        zLand = 0.0
        
        return Vector(xLand, yLand, zLand)
    else:
        return lastPosition
def getMaxHorizontalVel(self)
Expand source code
def getMaxHorizontalVel(self):
    maxHorizontalVel = 0
    for state in self.rigidBodyStates:
        horizVel = sqrt(state.velocity.X**2 + state.velocity.Y**2)
        if horizVel > maxHorizontalVel:
            maxHorizontalVel = horizVel
    return maxHorizontalVel
def getMaxSpeed(self)
Expand source code
def getMaxSpeed(self):
    maxVel = 0
    for state in self.rigidBodyStates:
        if state.velocity.length() > maxVel:
            maxVel = state.velocity.length()
    return maxVel
class SimDefinition (fileName=None, dictionary=None, disableDistributionSampling=False, silent=False, defaultDict=None, simDefParseStack=None)

Parse simulation definition files into a dictionary of string values accessible by string keys.

Inputs

  • fileName: (str) path to simulation definition file
  • dictionary: (dict[str,str]) Provide a pre-parsed dictionary equivalent to a simulation definition file - OVERRIDES fileName

  • disableDistributionSampling: (bool) Turn Monte Carlo sampling of normally-distributed parameters on/off

  • silent: (bool) Console output control
  • defaultDict: (dict[str,str] provide a custom dictionary of default values. If none is provided, defaultConfigValues is used.)
  • simDefParseStack: set(str) list of sim definition files in the current parse stack. Will throw an error if any of these files need to be loaded to generate the current sim definition # !include and !create [] from [] statements must form an acyclic graph of files to load (no circular loads)

Example

The file contents:
'SimControl{
    timeDiscretization RK4;
}'
Would be parsed into a single-key Python dictionary, stored in self.dict:
{ "SimControl.timeDiscretization": "RK4"}

Expand source code
class SimDefinition():

    #### Parsing / Initialization ####
    def __init__(self, fileName=None, dictionary=None, disableDistributionSampling=False, silent=False, defaultDict=None, simDefParseStack=None):
        '''
        Parse simulation definition files into a dictionary of string values accessible by string keys.

        Inputs:
            * fileName: (str) path to simulation definition file
            * dictionary: (dict[str,str]) Provide a pre-parsed dictionary equivalent to a simulation definition file - OVERRIDES fileName
            
            * disableDistributionSampling: (bool) Turn Monte Carlo sampling of normally-distributed parameters on/off
            * silent: (bool) Console output control
            * defaultDict: (dict[str,str] provide a custom dictionary of default values. If none is provided, defaultConfigValues is used.)
            * simDefParseStack: set(str) list of sim definition files in the current parse stack. Will throw an error if any of these files need to be loaded to generate the current sim definition
                # !include and !create [] from [] statements must form an acyclic graph of files to load (no circular loads)
        
        Example:
            The file contents:  
                'SimControl{  
                        timeDiscretization RK4;  
                }'  
            Would be parsed into a single-key Python dictionary, stored in self.dict:  
            `{ "SimControl.timeDiscretization": "RK4"}`
        
        '''
        self.silent = silent
        ''' Boolean, controls console output '''

        self.fileName = fileName

        self.disableDistributionSampling = disableDistributionSampling
        ''' Boolean - controls whether parameters which have standard deviations specified are actually sampled from a normal distribution. If True, the mean value is always returned. Chief use case for disabling sampling: Checking simulation convergence as the time step / target error is decreased. '''

        self.defaultDict = defaultConfigValues if (defaultDict == None) else defaultDict
        ''' Holds all of the defined default values. These will fill in for missing values in self.dict. Unless a different dictionary is specified, will hold a reference to `defaultConfigValues` '''
        
        self.monteCarloLogger = None 
        ''' Filled in by  Main.runMonteCarloSimulation() if running  Monte Carlo simulation. Type: `MAPLEAF.IO.Logging.MonteCarloLogger` '''

        self.dict = None # type: Dict[str:str]
        ''' Main dictionary of values, usually populated from a simulation definition file '''

        self.simDefParseStack = { self.fileName } if (simDefParseStack == None) else simDefParseStack
        ''' Keeps track of which files have already been loaded in the current parse stack. If these are loaded again we're in a parsing loop '''

        # Parse/Assign main values dictionary
        if dictionary != None:
            self.dict = dictionary
        elif fileName != None:
            self.dict = self._parseSimDefinitionFile(fileName)
        else:
            raise ValueError("No fileName or dictionary provided to initialize the SimDefinition")

        # Initialize tracking of default values used and unaccessed keys
        self._resetUsageTrackers()

        # Check if any probabilistic keys exist
        containsProbabilisticValues = False
        for key in self.dict:
            if "_stdDev" == key[-7:]:
                containsProbabilisticValues = True

        # Initialize instance of random.Random for Monte Carlo sampling
        if not disableDistributionSampling:
            try:
                randomSeed = self.getValue("MonteCarlo.randomSeed")
            except KeyError:
                randomSeed = random.randrange(1000000)
            
            if not silent and containsProbabilisticValues:
                print("Monte Carlo random seed: {}".format(randomSeed))
                
            self.rng = random.Random(randomSeed)
            ''' Instace of random.Random owned by this instance of SimDefinition. Random seed can be specified by the MonteCarlo.randomSeed parameter. Used for sampling all normal distributions for parameters that have std dev specified. '''

            self.resampleProbabilisticValues()            

    def _loadSubSimDefinition(self, path: str):
        ''' 
            In the parsing process, may need to load other sim definition files, use this function when doing that to detect circular references 
            path can be relative to the location of the current file, absolute, or relative to the MAPLEAF install directory
            
            Throws ValueError if circular parsing detected.
            Returns a new SimDefinition object
        '''
        filePath = getAbsoluteFilePath(path, str( Path(self.fileName).parent ))

        if filePath not in self.simDefParseStack:
            self.simDefParseStack.add(filePath)
            subSimDef = SimDefinition(filePath, simDefParseStack=self.simDefParseStack)
            self.simDefParseStack.remove(filePath)
            return subSimDef
        
        else:
            raise ValueError("Encountered circular reference while trying to parse SimDefinition file: {}, which references: {}, which is already in the current parse stack: {}".format(self.fileName, filePath, self.simDefParseStack))

    def _parseDictionaryContents(self, Dict, workingText, startLine: int, currDictName: str, allowKeyOverwriting=False) -> int:
        ''' 
            Parses an individual subdictionary in a simdefinition file.
            Calls itself recursively to parse further sub dictionaries.
            Saves parsed key-value pairs to Dict

            workingText should be of type list[str]

            Returns index of next line to parse
        '''
        i = startLine

        while i < len(workingText):
            line = workingText[i].strip()
            splitLine = line.split()
            
            if splitLine[0] == "!create":
                # Parse derived subdictionary
                i = self._parseDerivedDictionary(Dict, workingText, i, currDictName)

            elif splitLine[0] == "!include":
                # Include contents of another sim definition file
                filePath = line[line.index(" "):].strip() # Handle file names with spaces
                subDef = self._loadSubSimDefinition(filePath)

                # Add keys to current sim definition, inside current dictionary
                for subDefkey in subDef.dict:
                    if currDictName == "":
                        key = subDefkey
                    else:
                        key = currDictName + "." + subDefkey

                    Dict[key] = subDef.dict[subDefkey]

            elif line[-1] == '{':
                # Parse regular Subdictionary
                subDictName = line[:-1] # Remove whitespace and dict start bracket
                
                # Recursive call to parse subdictionary
                if currDictName == "":
                    i = self._parseDictionaryContents(Dict, workingText, i+1, subDictName, allowKeyOverwriting)
                else:
                    i = self._parseDictionaryContents(Dict, workingText, i+1, currDictName + "." + subDictName, allowKeyOverwriting)

            elif line == '}':
                #End current dictionary - continue parsing at next line
                return i
                        
            elif len(splitLine) > 1:
                # Save a space-separated key-value pair
                key = splitLine[0]
                value = " ".join(splitLine[1:])
                if currDictName == "":
                    keyString = key
                else:
                    keyString = currDictName + "." + key

                if not keyString in Dict or allowKeyOverwriting:
                    Dict[keyString] = value
                else:
                    raise ValueError("Duplicate Key: " + keyString + " in File: " + self.fileName)
            
            else:
                # Error: Line not recognized as a dict start/end or a key/value pair
                print(simDefinitionHelpMessage)
                raise ValueError("Problem parsing line {}: {}".format(i, line))

            # Next line
            i += 1

    def _parseDerivedDictionary(self, Dict, workingText, initializationLine: int, currDictName: str) -> int:
        '''
            Parse a 'derived' subdictionary, defined with the !create command in .mapleaf files

            Inputs:
                workingText: (list[str]) lines of text in .mapleaf file
                initializationLine: (int) index of line defining the derived dictionary to be parsed in workingText
                currDictName: (str) name of the dictionary containing the derived dictionary to be parsed. "" if at root level

            Returns:
                (int): index of the last line in the derived subdictionary
        '''
        # workingText[initializationLine] should be something like:
            # '    !create SubDictionary2 from Dictionary1.SubDictionary1{'
        definitionLine = shlex.split(workingText[initializationLine])

        # Figure out complete name of new dictionary
        if currDictName == '':
            derivedDictName = definitionLine[1]
        else:
            derivedDictName = currDictName + '.' + definitionLine[1]

        # Parent dict is last command. Remove opening curly bracket (last character), if present
        dictPath = definitionLine[-1][:-1] if ('{' == definitionLine[-1][-1]) else definitionLine[-1]

        #### Load Parent/Source (Sub)Dictionary ####
        if ":" in dictPath:
            # Importing dictionary from another file
            fileName = dictPath.split(":")[0]              
            subSimDef = self._loadSubSimDefinition(fileName)
            sourceDict = subSimDef.dict
            
            dictPath = dictPath.split(":")[1]
            keysInParentDict = subSimDef.getSubKeys(dictPath)                
        
        else:
            # Deriving from dictionary in current file
            # Get keys from parent dict
            keysInParentDict = self.getSubKeys(dictPath, Dict)
            sourceDict = Dict
        
        if len(keysInParentDict) == 0:
            raise ValueError("ERROR: Dictionary to derive from: {} is not defined before {} in {}.".format(dictPath, derivedDictName, self.fileName))
        
        # Fill out temporary dict, after applying all modifiers, add values to main Dict
        derivedDict = {}

        # Rename all the keys in the parentDict -> relocate them to the new (sub)dictionary
        for parentKey in keysInParentDict:
            key = parentKey.replace(dictPath, derivedDictName)
            derivedDict[key] = sourceDict[parentKey]


        #### Apply additional commands ####
        i = initializationLine + 1
        while i < len(workingText):
            line = workingText[i]
            command = shlex.split(line)

            def removeQuotes(string):
                string = string.replace("'", "")
                return string.replace('"', "")

            if command[0] == "!replace":
                # Get string to replace and its replacement
                toReplace = removeQuotes(command[1])
                replaceWith = removeQuotes(command[-1])

                derivedDictAfterReplace = {}
                for key in derivedDict:
                    newKey = key.replace(toReplace, replaceWith)
                    # .pop() gets the old value and also deletes it from the dictionary
                    newValue = derivedDict[key].replace(toReplace, replaceWith)
                    derivedDictAfterReplace[newKey] = newValue

                derivedDict = derivedDictAfterReplace

            elif command[0] == "!removeKeysContaining":
                stringToDelete = command[1]

                # Search for and remove any keys that contain stringToDelete
                keysToDelete = []
                for key in derivedDict:
                    if stringToDelete in key:
                        keysToDelete.append(key)

                for key in keysToDelete:
                    del derivedDict[key]                

            elif line[0] != "!":
                break # Done special commands - let the regular parser handle the rest

            else:
                raise ValueError("Command: {} not implemented. Try using !replace or !removeKeysContaining".format(line.split()[0]))

            i += 1

        #### Add derivedDict values to Dict ####
        for key in derivedDict:
            # Make sure we don't clobber existing values with poorly thought-out replace commands
            if key not in Dict:
                Dict[key] = derivedDict[key]
            else:
                raise ValueError("Derived dict key {} already exists".format(key, self.fileName))

        #### Parse any regular values in derived dict ####
        return self._parseDictionaryContents(Dict, workingText, i, derivedDictName, allowKeyOverwriting=True)

    def _replaceRelativeFilePathsWithAbsolutePaths(self, Dict):
        ''' 
            Tries to detect paths relative to the MAPLEAF installation directory and replaces them with absolute paths.
            This allows MAPLEAF to work when installed from pip and being run outside its installation directory.
        '''
        if self.fileName != None:
            fileDirectory = os.path.dirname(os.path.realpath(self.fileName))
        else:
            fileDirectory = None

        for key in Dict:
            # Iterate over all keys, looking for file path relative to the MAPLEAF repo
            val = Dict[key]
            
            # Remove leading dot/slash
            if val[:2] == "./":
                val = val[2:]

            if pathIsRelativeToRepository(val):
                # Replace the relative path with an absolute one
                Dict[key] = getAbsoluteFilePath(val)
            
            if isFileName(val):
                # Check if the file path is relative to the location of the simulation definition file
                if fileDirectory != None:
                    possibleLocation = os.path.join(fileDirectory, val)
                    if os.path.exists(possibleLocation):
                        Dict[key] = possibleLocation

    def _parseSimDefinitionFile(self, fileName):
        Dict = {}
        
        # Read all of the file's contents
        file = open(fileName, "r+")
        workingText = file.read()
        file.close()
        
        # Remove comments
        comment = re.compile("(?<!\\\)#.*") 
        workingText = re.sub(comment, "", workingText) 
        
        # Remove comment escape characters
        workingText = re.sub(r"\\(?=#)", "", workingText) 
        
        # Remove blank lines
        workingText = [line for line in workingText.split('\n') if line.strip() != '']
        
        # Start recursive parse by asking to parse the root-level dictionary
        self._parseDictionaryContents(Dict, workingText, 0, "")

        # Look for file paths relative to the MAPLEAF install location, replace them with absolute paths
        self._replaceRelativeFilePathsWithAbsolutePaths(Dict)

        return Dict

    #### Normal Usage ####
    def resampleProbabilisticValues(self, Dict=None):
        '''
            Normal Distribution Sampling:
                If (key + "_stdDev") exists and the value being returned is a scalar or Vector value, returns a scalar or vector sampled from a normal distribution
                    Where the mean of the normal distribution is taken to be the (original) value of 'key' (moved to 'key_mean' when this function first runs) and the standard deviation of the distribution is the value of 'key_stdDev'
                    For a vector value, a vector of standard deviations is expected
                For repeatable sampling, set the value "MonteCarlo.randomSeed" in the file loaded by this class
        '''
        if Dict is None:
            Dict = self.dict

        if not self.disableDistributionSampling:
            keys = list(Dict.keys()) # Get a list of keys at the beginning to avoid issues from the number of keys changing during iterations
            
            for key in keys:
                ### Sample any probabilistic values from normal distribution ###
                stdDevKey = key + "_stdDev"
                
                if stdDevKey in Dict:
                    logLine = None
                    meanKey = key + "_mean"

                    try:
                        meanString = Dict[meanKey]
                    except KeyError:
                        # Take the value of the variable as the mean if a _mean value is not provided
                        meanString = Dict[key]
                        Dict[meanKey] =  meanString

                    # Try parsing scalar values
                    try:
                        mu = float(meanString)
                        sigma = float(Dict[stdDevKey])

                        sampledValue = self.rng.gauss(mu, sigma)
                        Dict[key] = str(sampledValue)

                        logLine = "Sampling scalar parameter: {}, value: {:1.3f}".format(key, sampledValue)

                    except ValueError:
                        # Try parsing vector value
                        try:
                            muVec = Vector(meanString)
                            sigmaVec = Vector(Dict[stdDevKey])

                            sampledVec = Vector(*[ self.rng.gauss(mu, sigma) for mu, sigma in zip(muVec, sigmaVec)])
                            Dict[key] =  str(sampledVec)

                            logLine = "Sampling vector parameter: {}, value: ({:1.3f})".format(key, sampledVec)
                            
                        except ValueError:
                            # ValueError throws if either conversion to Vector fails
                            # Note that monte carlo / probabilistic variables can only be scalars or vectors
                            print("ERROR: Unable to parse probabilistic value: {} for key {} (or {} for key {}). Note that probabilistic values must be either scalars or vectors of length 3.".format(meanString, meanKey, self.getValue(stdDevKey), stdDevKey))
                            raise
                    
                    ### Logging ###
                    if logLine != None:
                        if self.monteCarloLogger != None:
                            self.monteCarloLogger.log(logLine)
                        elif not self.silent:
                            print(logLine)

    def getValue(self, key: str) -> str:
        """
            Input:
                Key should be a string of format "DictionaryName.SubdictionaryName.Key"
            Output:
                Always returns a string value
                Returns value from defaultConfigValues if key not present in current SimDefinition's dictionary
        """
        # Remove any whitespace from the key
        key = key.strip()

        ### Find string/mean value ###
        if self.dict.__contains__(key):
            if key in self.unaccessedFields: # Track which keys are accessed
                self.unaccessedFields.remove(key)
            return self.dict[key]
        
        elif key in self.defaultDict:
            self.defaultValuesUsed.add(key)
            return self.defaultDict[key]
        else:
            # Check if there's a class-based default value to return
            classBasedDefaultValue = self._getClassBasedDefaultValue(key)
            
            if classBasedDefaultValue != None:
                return classBasedDefaultValue
            else:
                raise KeyError("Key: " + key + " not found in {} or default config values".format(self.fileName))

    def setValue(self, key: str, value) -> None:
        '''
            Will add the entry if it's not present
        '''
        # The .strip() removes whitespace
        self.dict[key.strip()] = value

    def removeKey(self, key: str):
        if key in self.dict:
            return self.dict.pop(key)
        else:
            print("Warning: " + key + " not found, can't delete")
            return None

    def setIfAbsent(self, key: str, value):
        ''' Sets a value, only if it doesn't currently exist in the dictionary '''
        if not key in self.dict:
            self.setValue(key, value)

    def writeToFile(self, fileName: str, writeHeader=True) -> None:
        ''' 
            Write a (potentially modified) sim definition to file.
            Newly written file will not contain any comments! 
        '''
        self.fileName = fileName

        with open(fileName, 'w') as file:
            # Extract the fileName from the fileName variable, which may contain other folder names
            dictName = re.sub("^.*/", "", fileName)

            # Write Header
            if writeHeader:
                file.write("# MAPLEAF\n")
                file.write("# File: {}\n".format(fileName))
                file.write("# Autowritten on: " + str(datetime.now()) + "\n")

            # Sorting the keys before iterating through them ensures that dictionaries will be stored together
            sortedDict = sorted(self.dict.items())
            currDicts = []
            for key in sortedDict:
                key = key[0]
                dicts = key.split('.')[:-1]

                # Need to get be in the appropriate dictionary before writing the key, value pair
                if dicts != currDicts:
                    
                    #Close any uneeded dictionaries
                    dictDepth = currDicts.__len__()
                    while dictDepth > 0:
                        if dictDepth > dicts.__len__():
                            file.write("\t"*(dictDepth-1) + "}\n")
                        elif currDicts[dictDepth-1] != dicts[dictDepth-1]:
                            file.write("\t"*(dictDepth-1) + "}\n")
                        else:
                            break
                        
                        dictDepth = dictDepth - 1

                    openedNewDict = False

                    #Open any new dictionaries
                    while dictDepth < dicts.__len__():
                        newDict = dicts[dictDepth]
                        file.write("\n" + "\t" * dictDepth + newDict + "{\n")
                        dictDepth = dictDepth + 1
                        openedNewDict = True
                    
                    if not openedNewDict:
                        # If no new dictionary was openend after closing unneeded ones, add a spacing line before writing keys/values
                        file.write("\n")

                    currDicts = dicts

                #Add the key, value
                dictDepth = currDicts.__len__()
                realKey = re.sub("^([^\.]*\.)+", "", key)
                file.write( "\t"*dictDepth + realKey + "\t" + self.dict[key] + "\n")

            #Close any open dictionaries
            dictDepth = currDicts.__len__()
            while dictDepth > 0:
                dictDepth = dictDepth - 1
                file.write("\t"*dictDepth + "}\n")

    #### Introspection / Key Gymnastics ####
    def findKeysContaining(self, keyContains: List[str]) -> List[str]:
        '''
            Returns a list of all keys that contain any of the strings in keyContains
            
            ## Example  
                findKeysContaining(["class"]) ->  
                [ "Rocket.class", "Rocket.Sustainer.class", "Rocket.Sustainer.Nosecone.class", etc... ]
        '''
        if not isinstance(keyContains, list):
            keyContains = [ keyContains ]
        
        matchingKeys = []
        for key in self.dict.keys():
            match = True
            for str in keyContains:
                if str not in key:
                    match = False
                    break
            
            if match:
                matchingKeys.append(key)
        
        if len(matchingKeys) > 0:
            return matchingKeys
        else:
            return None

    def getSubKeys(self, key: str, Dict=None) -> List[str]:
        '''
            Returns a list of all keys that are children of key

            ## Example  
                getSubKeys("Rocket") ->  
                [ "Rocket.position", "Rocket.Sustainer.NoseCone.mass", "Rocket.Sustainer.RecoverySystem.position", etc... ]
        '''
        #TODO: Improve speed by keeping dict sorted, then use binary search to locate first/last subkeys
        Dict = self.dict if (Dict == None) else Dict

        subKeys = []
        for currentKey in Dict.keys():
            if isSubKey(key, currentKey):
                subKeys.append(currentKey)
        
        return subKeys

    def getImmediateSubKeys(self, key: str) -> List[str]:
        """ 
            Returns all keys that are immediate children of the parentKey (one 'level' lower)
            
            .. note:: Will not return subdictionaries, only keys that have a value associated with them. Use self.getImmediateSubDicts() to discover sub-dictionaries

            ## Example:
                getImmediateSubKeys("Rocket") ->  
                [ "Rocket.name", "Rocket.position", "Rocket.velocity", etc...]
        """
        results = set()
        for potentialChildKey in self.dict.keys():
            # Iterate through all keys - check if they are children of currentPath
            if isSubKey(key, potentialChildKey):
                # If so, get the part of the key that is the immediate child of currentPath
                immediateSubkey = getImmediateSubKey(key, potentialChildKey)
                
                if immediateSubkey in self.dict:
                    # If we haven't got it already, save it
                    results.add(immediateSubkey)

        return list(results)

    def getImmediateSubDicts(self, key: str) -> List[str]:
        '''
            Returns list of names of immediate subdictionaries

            ## Example
                getImmediateSubDicts("Rocket") ->
                [ "Rocket.StageOne", "Rocket.StageTwo", "Rocket.ControlSystem", etc... ]

            .. note:: This example would not return a dictionry like: "Rocket.StageOne.FinSet" because it's not an immediate subdictionary of "Rocket"
        '''
        keyLevel = getKeyLevel(key)
        subKeys = self.getSubKeys(key)

        subDictionaries = set()
        for subKey in subKeys:
            subKeyLevel = getKeyLevel(subKey)
            if subKeyLevel - keyLevel > 1:
                # A subkey would have 1 level higher
                # A subkey of a subdictionary would have 2 levels higher - this is what we're looking for
                subDictKey = getParentKeyAtLevel(subKey, keyLevel+1)
                subDictionaries.add(subDictKey)
        
        return list(subDictionaries)

    def _getClassBasedDefaultValue(self, key: str) -> Union[str, None]:
        ''' 
            Returns class-based default value from defaultConfigValues if it exists. Otherwise returns None 
            
            Will attempt to find class-based default values for every longer prefixes of a key:
                key = "Rocket.Sustainer.canards.trailingEdge.shape"
                Attempt1 = "Rocket.Sustainer.canards.trailingEdge.class" -> Fail
                Attempt2 = "Rocket.Sustainer.canards.class" -> FinSet -> look up 'FinSet.trailingEdge.shape' in defaultDict -> if there, return it, otherwise return None
        '''
        splitLevel = getKeyLevel(key)

        while splitLevel >= 0:
            prefix, suffix = splitKeyAtLevel(key, splitLevel)
            
            try:
                classKey = prefix + ".class"
                className = self.dict[classKey]                

                # As soon as we arrive at an item with a class, search terminates
                try:
                    classBasedDefaultKey = className + "." + suffix
                    defaultValue = self.defaultDict[classBasedDefaultKey]

                    # Track that we've used a default value
                    self.defaultValuesUsed.add(classBasedDefaultKey)
                    
                    # if the classKey was useful, count it as 'used'
                    if classKey in self.unaccessedFields: 
                        self.unaccessedFields.remove(classKey)
                        
                    return defaultValue
                except KeyError:
                    return None # class-based default value not found
            
            except KeyError:
                pass # prefix.class not present

            # Move one level up the dictionary for next attempt
            splitLevel -= 1
        
        return None

    #### Usage Reporting ####
    def printUnusedKeys(self):
        '''
            Checks which keys in the present simulation definition have not yet been accessed.
            Prints a list of those to the console.
        '''
        if len(self.unaccessedFields) > 0:
            print("\nWarning: The following keys were loaded from: {} but never accessed:".format(self.fileName))
            for key in sorted(self.unaccessedFields):
                value = self.dict[key]
                print("{:<45}{}".format(key+":", value))
            print("")

    def printDefaultValuesUsed(self):
        '''
            Checks which default values have been used since the creation of the current instance of SimDefinition. Prints those to the console.
        '''
        if len(self.defaultValuesUsed):
            print("\nWarning: The following default values were used in this simulation:")
            for key in sorted(self.defaultValuesUsed):
                value = self.defaultDict[key]
                print("{:<45}{}".format(key+":", value))
            print("\nIf this was not intended, override the default values by adding the above information to your simulation definition file.\n")
        
    def _resetUsageTrackers(self):
        # Create a dictionary to keep track of which attributed have been accessed (initially none)
        self.unaccessedFields = set(self.dict.keys())
        # Create a list to track which default values have been used
        self.defaultValuesUsed = set()

    #### Utilities ####
    def __str__(self):
        result = ""
        result += "File: " + self.fileName + "\n"

        for key, value in self.dict.items():
            result += "{}: {}\n".format(key, value)

        result += "\n"

        return result

    def __eq__ (self, simDef2):
        try:
            if self.dict == simDef2.dict:
                return True
            else:
                return False
        except AttributeError:
            return False

    def __contains__(self, key):
        ''' Only checks whether 'key' was parsed from the file. Ignores default values '''
        return key in self.dict

Instance variables

var defaultDict

Holds all of the defined default values. These will fill in for missing values in self.dict. Unless a different dictionary is specified, will hold a reference to defaultConfigValues

var dict

Main dictionary of values, usually populated from a simulation definition file

var disableDistributionSampling

Boolean - controls whether parameters which have standard deviations specified are actually sampled from a normal distribution. If True, the mean value is always returned. Chief use case for disabling sampling: Checking simulation convergence as the time step / target error is decreased.

var monteCarloLogger

Filled in by Main.runMonteCarloSimulation() if running Monte Carlo simulation. Type: MonteCarloLogger

var silent

Boolean, controls console output

var simDefParseStack

Keeps track of which files have already been loaded in the current parse stack. If these are loaded again we're in a parsing loop

Methods

def findKeysContaining(self, keyContains: List[str]) ‑> List[str]

Returns a list of all keys that contain any of the strings in keyContains

Example

findKeysContaining(["class"]) ->  
[ "Rocket.class", "Rocket.Sustainer.class", "Rocket.Sustainer.Nosecone.class", etc... ]
Expand source code
def findKeysContaining(self, keyContains: List[str]) -> List[str]:
    '''
        Returns a list of all keys that contain any of the strings in keyContains
        
        ## Example  
            findKeysContaining(["class"]) ->  
            [ "Rocket.class", "Rocket.Sustainer.class", "Rocket.Sustainer.Nosecone.class", etc... ]
    '''
    if not isinstance(keyContains, list):
        keyContains = [ keyContains ]
    
    matchingKeys = []
    for key in self.dict.keys():
        match = True
        for str in keyContains:
            if str not in key:
                match = False
                break
        
        if match:
            matchingKeys.append(key)
    
    if len(matchingKeys) > 0:
        return matchingKeys
    else:
        return None
def getImmediateSubDicts(self, key: str) ‑> List[str]

Returns list of names of immediate subdictionaries

Example

getImmediateSubDicts("Rocket") ->
[ "Rocket.StageOne", "Rocket.StageTwo", "Rocket.ControlSystem", etc... ]

Note: This example would not return a dictionry like: "Rocket.StageOne.FinSet" because it's not an immediate subdictionary of "Rocket"

Expand source code
def getImmediateSubDicts(self, key: str) -> List[str]:
    '''
        Returns list of names of immediate subdictionaries

        ## Example
            getImmediateSubDicts("Rocket") ->
            [ "Rocket.StageOne", "Rocket.StageTwo", "Rocket.ControlSystem", etc... ]

        .. note:: This example would not return a dictionry like: "Rocket.StageOne.FinSet" because it's not an immediate subdictionary of "Rocket"
    '''
    keyLevel = getKeyLevel(key)
    subKeys = self.getSubKeys(key)

    subDictionaries = set()
    for subKey in subKeys:
        subKeyLevel = getKeyLevel(subKey)
        if subKeyLevel - keyLevel > 1:
            # A subkey would have 1 level higher
            # A subkey of a subdictionary would have 2 levels higher - this is what we're looking for
            subDictKey = getParentKeyAtLevel(subKey, keyLevel+1)
            subDictionaries.add(subDictKey)
    
    return list(subDictionaries)
def getImmediateSubKeys(self, key: str) ‑> List[str]

Returns all keys that are immediate children of the parentKey (one 'level' lower)

Note: Will not return subdictionaries, only keys that have a value associated with them. Use self.getImmediateSubDicts() to discover sub-dictionaries

Example:

getImmediateSubKeys("Rocket") ->  
[ "Rocket.name", "Rocket.position", "Rocket.velocity", etc...]
Expand source code
def getImmediateSubKeys(self, key: str) -> List[str]:
    """ 
        Returns all keys that are immediate children of the parentKey (one 'level' lower)
        
        .. note:: Will not return subdictionaries, only keys that have a value associated with them. Use self.getImmediateSubDicts() to discover sub-dictionaries

        ## Example:
            getImmediateSubKeys("Rocket") ->  
            [ "Rocket.name", "Rocket.position", "Rocket.velocity", etc...]
    """
    results = set()
    for potentialChildKey in self.dict.keys():
        # Iterate through all keys - check if they are children of currentPath
        if isSubKey(key, potentialChildKey):
            # If so, get the part of the key that is the immediate child of currentPath
            immediateSubkey = getImmediateSubKey(key, potentialChildKey)
            
            if immediateSubkey in self.dict:
                # If we haven't got it already, save it
                results.add(immediateSubkey)

    return list(results)
def getSubKeys(self, key: str, Dict=None) ‑> List[str]

Returns a list of all keys that are children of key

Example

getSubKeys("Rocket") ->  
[ "Rocket.position", "Rocket.Sustainer.NoseCone.mass", "Rocket.Sustainer.RecoverySystem.position", etc... ]
Expand source code
def getSubKeys(self, key: str, Dict=None) -> List[str]:
    '''
        Returns a list of all keys that are children of key

        ## Example  
            getSubKeys("Rocket") ->  
            [ "Rocket.position", "Rocket.Sustainer.NoseCone.mass", "Rocket.Sustainer.RecoverySystem.position", etc... ]
    '''
    #TODO: Improve speed by keeping dict sorted, then use binary search to locate first/last subkeys
    Dict = self.dict if (Dict == None) else Dict

    subKeys = []
    for currentKey in Dict.keys():
        if isSubKey(key, currentKey):
            subKeys.append(currentKey)
    
    return subKeys
def getValue(self, key: str) ‑> str

Input

Key should be a string of format "DictionaryName.SubdictionaryName.Key"

Output

Always returns a string value Returns value from defaultConfigValues if key not present in current SimDefinition's dictionary

Expand source code
def getValue(self, key: str) -> str:
    """
        Input:
            Key should be a string of format "DictionaryName.SubdictionaryName.Key"
        Output:
            Always returns a string value
            Returns value from defaultConfigValues if key not present in current SimDefinition's dictionary
    """
    # Remove any whitespace from the key
    key = key.strip()

    ### Find string/mean value ###
    if self.dict.__contains__(key):
        if key in self.unaccessedFields: # Track which keys are accessed
            self.unaccessedFields.remove(key)
        return self.dict[key]
    
    elif key in self.defaultDict:
        self.defaultValuesUsed.add(key)
        return self.defaultDict[key]
    else:
        # Check if there's a class-based default value to return
        classBasedDefaultValue = self._getClassBasedDefaultValue(key)
        
        if classBasedDefaultValue != None:
            return classBasedDefaultValue
        else:
            raise KeyError("Key: " + key + " not found in {} or default config values".format(self.fileName))
def printDefaultValuesUsed(self)

Checks which default values have been used since the creation of the current instance of SimDefinition. Prints those to the console.

Expand source code
def printDefaultValuesUsed(self):
    '''
        Checks which default values have been used since the creation of the current instance of SimDefinition. Prints those to the console.
    '''
    if len(self.defaultValuesUsed):
        print("\nWarning: The following default values were used in this simulation:")
        for key in sorted(self.defaultValuesUsed):
            value = self.defaultDict[key]
            print("{:<45}{}".format(key+":", value))
        print("\nIf this was not intended, override the default values by adding the above information to your simulation definition file.\n")
def printUnusedKeys(self)

Checks which keys in the present simulation definition have not yet been accessed. Prints a list of those to the console.

Expand source code
def printUnusedKeys(self):
    '''
        Checks which keys in the present simulation definition have not yet been accessed.
        Prints a list of those to the console.
    '''
    if len(self.unaccessedFields) > 0:
        print("\nWarning: The following keys were loaded from: {} but never accessed:".format(self.fileName))
        for key in sorted(self.unaccessedFields):
            value = self.dict[key]
            print("{:<45}{}".format(key+":", value))
        print("")
def removeKey(self, key: str)
Expand source code
def removeKey(self, key: str):
    if key in self.dict:
        return self.dict.pop(key)
    else:
        print("Warning: " + key + " not found, can't delete")
        return None
def resampleProbabilisticValues(self, Dict=None)

Normal Distribution Sampling: If (key + "_stdDev") exists and the value being returned is a scalar or Vector value, returns a scalar or vector sampled from a normal distribution Where the mean of the normal distribution is taken to be the (original) value of 'key' (moved to 'key_mean' when this function first runs) and the standard deviation of the distribution is the value of 'key_stdDev' For a vector value, a vector of standard deviations is expected For repeatable sampling, set the value "MonteCarlo.randomSeed" in the file loaded by this class

Expand source code
def resampleProbabilisticValues(self, Dict=None):
    '''
        Normal Distribution Sampling:
            If (key + "_stdDev") exists and the value being returned is a scalar or Vector value, returns a scalar or vector sampled from a normal distribution
                Where the mean of the normal distribution is taken to be the (original) value of 'key' (moved to 'key_mean' when this function first runs) and the standard deviation of the distribution is the value of 'key_stdDev'
                For a vector value, a vector of standard deviations is expected
            For repeatable sampling, set the value "MonteCarlo.randomSeed" in the file loaded by this class
    '''
    if Dict is None:
        Dict = self.dict

    if not self.disableDistributionSampling:
        keys = list(Dict.keys()) # Get a list of keys at the beginning to avoid issues from the number of keys changing during iterations
        
        for key in keys:
            ### Sample any probabilistic values from normal distribution ###
            stdDevKey = key + "_stdDev"
            
            if stdDevKey in Dict:
                logLine = None
                meanKey = key + "_mean"

                try:
                    meanString = Dict[meanKey]
                except KeyError:
                    # Take the value of the variable as the mean if a _mean value is not provided
                    meanString = Dict[key]
                    Dict[meanKey] =  meanString

                # Try parsing scalar values
                try:
                    mu = float(meanString)
                    sigma = float(Dict[stdDevKey])

                    sampledValue = self.rng.gauss(mu, sigma)
                    Dict[key] = str(sampledValue)

                    logLine = "Sampling scalar parameter: {}, value: {:1.3f}".format(key, sampledValue)

                except ValueError:
                    # Try parsing vector value
                    try:
                        muVec = Vector(meanString)
                        sigmaVec = Vector(Dict[stdDevKey])

                        sampledVec = Vector(*[ self.rng.gauss(mu, sigma) for mu, sigma in zip(muVec, sigmaVec)])
                        Dict[key] =  str(sampledVec)

                        logLine = "Sampling vector parameter: {}, value: ({:1.3f})".format(key, sampledVec)
                        
                    except ValueError:
                        # ValueError throws if either conversion to Vector fails
                        # Note that monte carlo / probabilistic variables can only be scalars or vectors
                        print("ERROR: Unable to parse probabilistic value: {} for key {} (or {} for key {}). Note that probabilistic values must be either scalars or vectors of length 3.".format(meanString, meanKey, self.getValue(stdDevKey), stdDevKey))
                        raise
                
                ### Logging ###
                if logLine != None:
                    if self.monteCarloLogger != None:
                        self.monteCarloLogger.log(logLine)
                    elif not self.silent:
                        print(logLine)
def setIfAbsent(self, key: str, value)

Sets a value, only if it doesn't currently exist in the dictionary

Expand source code
def setIfAbsent(self, key: str, value):
    ''' Sets a value, only if it doesn't currently exist in the dictionary '''
    if not key in self.dict:
        self.setValue(key, value)
def setValue(self, key: str, value) ‑> None

Will add the entry if it's not present

Expand source code
def setValue(self, key: str, value) -> None:
    '''
        Will add the entry if it's not present
    '''
    # The .strip() removes whitespace
    self.dict[key.strip()] = value
def writeToFile(self, fileName: str, writeHeader=True) ‑> None

Write a (potentially modified) sim definition to file. Newly written file will not contain any comments!

Expand source code
def writeToFile(self, fileName: str, writeHeader=True) -> None:
    ''' 
        Write a (potentially modified) sim definition to file.
        Newly written file will not contain any comments! 
    '''
    self.fileName = fileName

    with open(fileName, 'w') as file:
        # Extract the fileName from the fileName variable, which may contain other folder names
        dictName = re.sub("^.*/", "", fileName)

        # Write Header
        if writeHeader:
            file.write("# MAPLEAF\n")
            file.write("# File: {}\n".format(fileName))
            file.write("# Autowritten on: " + str(datetime.now()) + "\n")

        # Sorting the keys before iterating through them ensures that dictionaries will be stored together
        sortedDict = sorted(self.dict.items())
        currDicts = []
        for key in sortedDict:
            key = key[0]
            dicts = key.split('.')[:-1]

            # Need to get be in the appropriate dictionary before writing the key, value pair
            if dicts != currDicts:
                
                #Close any uneeded dictionaries
                dictDepth = currDicts.__len__()
                while dictDepth > 0:
                    if dictDepth > dicts.__len__():
                        file.write("\t"*(dictDepth-1) + "}\n")
                    elif currDicts[dictDepth-1] != dicts[dictDepth-1]:
                        file.write("\t"*(dictDepth-1) + "}\n")
                    else:
                        break
                    
                    dictDepth = dictDepth - 1

                openedNewDict = False

                #Open any new dictionaries
                while dictDepth < dicts.__len__():
                    newDict = dicts[dictDepth]
                    file.write("\n" + "\t" * dictDepth + newDict + "{\n")
                    dictDepth = dictDepth + 1
                    openedNewDict = True
                
                if not openedNewDict:
                    # If no new dictionary was openend after closing unneeded ones, add a spacing line before writing keys/values
                    file.write("\n")

                currDicts = dicts

            #Add the key, value
            dictDepth = currDicts.__len__()
            realKey = re.sub("^([^\.]*\.)+", "", key)
            file.write( "\t"*dictDepth + realKey + "\t" + self.dict[key] + "\n")

        #Close any open dictionaries
        dictDepth = currDicts.__len__()
        while dictDepth > 0:
            dictDepth = dictDepth - 1
            file.write("\t"*dictDepth + "}\n")
class SubDictReader (stringPathToThisItemsSubDictionary, simDefinition)

Example stringPathToThisItemsSubDictionary = 'Rocket.Stage1.Nosecone' if we're initializing a nosecone object

Expand source code
class SubDictReader():

    def __init__(self, stringPathToThisItemsSubDictionary, simDefinition):
        '''
            Example stringPathToThisItemsSubDictionary = 'Rocket.Stage1.Nosecone' if we're initializing a nosecone object
        '''
        self.simDefDictPathToReadFrom = stringPathToThisItemsSubDictionary
        self.simDefinition = simDefinition

    def getString(self, key):
        ''' 
            Pass in either relative key or absolute key:
                Ex 1 (Relative): If object subdictionary (self.simDefDictPathToReadFrom) is 'Rocket.Stage1.Nosecone', relative keys could be 'position' or 'aspectRatio'
                    These would retrieve Rocket.Stage1.Nosecone.position or Rocket.Stage1.Nosecone.aspectRatio from the sim definition
                Ex 2 (Absolute): Can also pass in full absolute key, like 'Rocket.Stage1.Nosecone.position', and it will retrieve that value, as long as there isn't a 'path collision' with a relative path
        '''
        try:
            return self.simDefinition.getValue(self.simDefDictPathToReadFrom + "." + key)
        except KeyError:
            try:
                return self.simDefinition.getValue(key)
            except KeyError:
                attemptedKey1 = self.simDefDictPathToReadFrom + "." + key
                attemptedKey2 = key
                raise KeyError("{} and {} not found in {} or in default value dictionary".format(attemptedKey1, attemptedKey2, self.simDefinition.fileName))

    #### Get parsed values ####
    def getInt(self, key: str) -> int:
        return int(self.getString(key))

    def getFloat(self, key: str) -> float:
        return float(self.getString(key))

    def getVector(self, key: str) -> Vector:
        return Vector(self.getString(key))

    def getBool(self, key: str) -> bool:
        return strtobool(self.getString(key))

    #### Try get values (return specified default value if not found) ####
    def tryGetString(self, key: str, defaultValue: Union[None, str]=None):
        try:
            return self.getString(key)
        except KeyError:
            return defaultValue

    def tryGetInt(self, key: str, defaultValue: Union[None, int]=None):
        try:
            return int(self.getString(key))
        except KeyError:
            return defaultValue

    def tryGetFloat(self, key: str, defaultValue: Union[None, float]=None):
        try:
            return float(self.getString(key))
        except KeyError:
            return defaultValue

    def tryGetVector(self, key: str, defaultValue: Union[None, Vector]=None):
        try:
            return Vector(self.getString(key))
        except KeyError:
            return defaultValue

    def tryGetBool(self, key: str, defaultValue: Union[None, bool]=None):
        try:
            return strtobool(self.getString(key))
        except KeyError:
            return defaultValue

    #### Introspection ####
    def getImmediateSubDicts(self, key=None) -> List[str]:
        if key == None:
            key = self.simDefDictPathToReadFrom
        else:
            absKey = self.simDefDictPathToReadFrom + "." + key
            result = self.simDefinition.getImmediateSubDicts(absKey)
        
            if len(result) != 0: # If no subdicts found, try other key
                return result
        
        return self.simDefinition.getImmediateSubDicts(key)

    def getSubKeys(self, key=None) -> List[str]:
        if key == None:
            key = self.simDefDictPathToReadFrom
        else:
            absKey = self.simDefDictPathToReadFrom + "." + key
            result = self.simDefinition.getSubKeys(absKey)

            if len(result) != 0: # If nothing found, try another key
                return result

        return self.simDefinition.getSubKeys(key)

    def getImmediateSubKeys(self, key=None) -> List[str]:
        if key == None:
            key = self.simDefDictPathToReadFrom
        else:
            absKey = self.simDefDictPathToReadFrom + "." + key
            result = self.simDefinition.getImmediateSubKeys(absKey)

            if len(result) != 0: # If nothing found, try another key
                return result

        return self.simDefinition.getImmediateSubKeys(key)

    def getDictName(self) -> str:
        lastDotIndex = self.simDefDictPathToReadFrom.rfind('.')
        return self.simDefDictPathToReadFrom[lastDotIndex+1:]

Subclasses

Methods

def getBool(self, key: str) ‑> bool
Expand source code
def getBool(self, key: str) -> bool:
    return strtobool(self.getString(key))
def getDictName(self) ‑> str
Expand source code
def getDictName(self) -> str:
    lastDotIndex = self.simDefDictPathToReadFrom.rfind('.')
    return self.simDefDictPathToReadFrom[lastDotIndex+1:]
def getFloat(self, key: str) ‑> float
Expand source code
def getFloat(self, key: str) -> float:
    return float(self.getString(key))
def getImmediateSubDicts(self, key=None) ‑> List[str]
Expand source code
def getImmediateSubDicts(self, key=None) -> List[str]:
    if key == None:
        key = self.simDefDictPathToReadFrom
    else:
        absKey = self.simDefDictPathToReadFrom + "." + key
        result = self.simDefinition.getImmediateSubDicts(absKey)
    
        if len(result) != 0: # If no subdicts found, try other key
            return result
    
    return self.simDefinition.getImmediateSubDicts(key)
def getImmediateSubKeys(self, key=None) ‑> List[str]
Expand source code
def getImmediateSubKeys(self, key=None) -> List[str]:
    if key == None:
        key = self.simDefDictPathToReadFrom
    else:
        absKey = self.simDefDictPathToReadFrom + "." + key
        result = self.simDefinition.getImmediateSubKeys(absKey)

        if len(result) != 0: # If nothing found, try another key
            return result

    return self.simDefinition.getImmediateSubKeys(key)
def getInt(self, key: str) ‑> int
Expand source code
def getInt(self, key: str) -> int:
    return int(self.getString(key))
def getString(self, key)

Pass in either relative key or absolute key: Ex 1 (Relative): If object subdictionary (self.simDefDictPathToReadFrom) is 'Rocket.Stage1.Nosecone', relative keys could be 'position' or 'aspectRatio' These would retrieve Rocket.Stage1.Nosecone.position or Rocket.Stage1.Nosecone.aspectRatio from the sim definition Ex 2 (Absolute): Can also pass in full absolute key, like 'Rocket.Stage1.Nosecone.position', and it will retrieve that value, as long as there isn't a 'path collision' with a relative path

Expand source code
def getString(self, key):
    ''' 
        Pass in either relative key or absolute key:
            Ex 1 (Relative): If object subdictionary (self.simDefDictPathToReadFrom) is 'Rocket.Stage1.Nosecone', relative keys could be 'position' or 'aspectRatio'
                These would retrieve Rocket.Stage1.Nosecone.position or Rocket.Stage1.Nosecone.aspectRatio from the sim definition
            Ex 2 (Absolute): Can also pass in full absolute key, like 'Rocket.Stage1.Nosecone.position', and it will retrieve that value, as long as there isn't a 'path collision' with a relative path
    '''
    try:
        return self.simDefinition.getValue(self.simDefDictPathToReadFrom + "." + key)
    except KeyError:
        try:
            return self.simDefinition.getValue(key)
        except KeyError:
            attemptedKey1 = self.simDefDictPathToReadFrom + "." + key
            attemptedKey2 = key
            raise KeyError("{} and {} not found in {} or in default value dictionary".format(attemptedKey1, attemptedKey2, self.simDefinition.fileName))
def getSubKeys(self, key=None) ‑> List[str]
Expand source code
def getSubKeys(self, key=None) -> List[str]:
    if key == None:
        key = self.simDefDictPathToReadFrom
    else:
        absKey = self.simDefDictPathToReadFrom + "." + key
        result = self.simDefinition.getSubKeys(absKey)

        if len(result) != 0: # If nothing found, try another key
            return result

    return self.simDefinition.getSubKeys(key)
def getVector(self, key: str) ‑> Vector
Expand source code
def getVector(self, key: str) -> Vector:
    return Vector(self.getString(key))
def tryGetBool(self, key: str, defaultValue: Optional[None] = None)
Expand source code
def tryGetBool(self, key: str, defaultValue: Union[None, bool]=None):
    try:
        return strtobool(self.getString(key))
    except KeyError:
        return defaultValue
def tryGetFloat(self, key: str, defaultValue: Optional[None] = None)
Expand source code
def tryGetFloat(self, key: str, defaultValue: Union[None, float]=None):
    try:
        return float(self.getString(key))
    except KeyError:
        return defaultValue
def tryGetInt(self, key: str, defaultValue: Optional[None] = None)
Expand source code
def tryGetInt(self, key: str, defaultValue: Union[None, int]=None):
    try:
        return int(self.getString(key))
    except KeyError:
        return defaultValue
def tryGetString(self, key: str, defaultValue: Optional[None] = None)
Expand source code
def tryGetString(self, key: str, defaultValue: Union[None, str]=None):
    try:
        return self.getString(key)
    except KeyError:
        return defaultValue
def tryGetVector(self, key: str, defaultValue: Optional[None] = None)
Expand source code
def tryGetVector(self, key: str, defaultValue: Union[None, Vector]=None):
    try:
        return Vector(self.getString(key))
    except KeyError:
        return defaultValue
class TimeStepLog (...)

Adds functionality for post processing specific to time step logs for mapleaf

columnNames should be a list of strings To get references to the log columns, initialize the log as empty and then use the addColumn(s) functions. Those return the list(s) they add to the log for direct access

Ancestors

Methods

def writeToCSV(...)

Inherited members