TurbulenceConvectionUtils
CalibrateEDMF.TurbulenceConvectionUtils.ModelEvaluator — TypeModelEvaluatorA structure containing the information required to perform a forward model evaluation at a given parameter vector.
Fields
- param_cons::Vector{FT} where FT<:Real: Parameter vector in constrained (i.e. physical) space
- param_names::Vector{String}: Parameter names associated with parameter vector
- param_map::CalibrateEDMF.HelperFuncs.ParameterMap: A mapping operator to define relations between parameters. See- ParameterMapfor details.
- ref_models::Vector{CalibrateEDMF.ReferenceModels.ReferenceModel}: Vector of reference models
- ref_stats::CalibrateEDMF.ReferenceStats.ReferenceStatistics: Reference statistics for the inverse problem
Constructors
ModelEvaluator(
;
    param_cons,
    param_names,
    param_map,
    ref_models,
    ref_stats
)defined at /home/runner/work/CalibrateEDMF.jl/CalibrateEDMF.jl/src/TurbulenceConvectionUtils.jl:51.
ModelEvaluator(param_cons, param_names, param_map, RM, RS)defined at /home/runner/work/CalibrateEDMF.jl/CalibrateEDMF.jl/src/TurbulenceConvectionUtils.jl:64.
CalibrateEDMF.TurbulenceConvectionUtils.run_SCM — Functionrun_SCM(
    u::Vector{FT},
    u_names::Vector{String},
    param_map::ParameterMap,
    RM::Vector{ReferenceModel},
    RS::ReferenceStatistics;
    error_check::Bool = false,
    failure_handler = "high_loss",
) where {FT <: Real}
run_SCM(
    ME::ModelEvaluator;
    error_check::Bool = false,
    failure_handler = "high_loss",
) where {FT <: Real}Run the single-column model (SCM) using a set of parameters and return the value of outputs defined in y_names, possibly after normalization and projection onto lower dimensional space using PCA.
The function also outputs a boolean diagnosing whether the requested SCM simulation failed.
Inputs:
- u:: Values of parameters to be used in simulations.
- u_names:: SCM names for parameters- u.
- param_map:: A mapping operator to define relations between parameters. See- ParameterMapfor details.
- RM:: Vector of- ReferenceModels
- RS:: reference statistics for simulation
- error_check:: Returns as an additional argument whether the SCM call errored.
- failure_handler:: Method used to handle failed simulations.
Outputs:
- sim_dirs:: Vector of simulation output directories
- g_scm:: Vector of model evaluations concatenated for all flow configurations.
- g_scm_pca:: Projection of- g_scmonto principal subspace spanned by eigenvectors.
- model_error:: Whether the simulation errored with the requested configuration.
CalibrateEDMF.TurbulenceConvectionUtils.run_reference_SCM — Functionrun_reference_SCM(m::ReferenceModel; overwrite::Bool = false, run_single_timestep = true)Run the single-column model (SCM) for a reference model object using default parameters, and write the output to file.
Inputs:
- m:: A- ReferenceModel.
- overwrite:: if true, run TC.jl and overwrite existing simulation files.
- run_single_timestep:: if true, run only one time step.
CalibrateEDMF.TurbulenceConvectionUtils.eval_single_ref_model — Functioneval_single_ref_model(
    m_index::IT,
    m::ReferenceModel,
    RS::ReferenceStatistics,
    u::Vector{FT},
    u_names::Vector{String},
    param_map::ParameterMap,
) where {FT <: Real, IT <: Int}Runs the single-column model (SCM) under a single configuration (i.e., ReferenceModel) using a set of parameters u, and returns the forward model evaluation in both the original and the latent PCA space.
Inputs:
- m_index:: The index of the ReferenceModel within the overarching ref_models vector used to construct the ReferenceStatistics.
- m:: A ReferenceModel.
- RS:: reference statistics for simulation
- u:: Values of parameters to be used in simulations.
- u_names:: SCM names for parameters- u.
- param_map:: A mapping operator to define relations between parameters. See- ParameterMapfor details.
Outputs:
- sim_dir:: Simulation output directory.
- g_scm:: Forward model evaluation in original output space.
- g_scm_pca:: Projection of- g_scmonto principal subspace spanned by eigenvectors.
- model_error:: Whether the simulation errored with the requested configuration.
CalibrateEDMF.TurbulenceConvectionUtils.run_SCM_handler — Functionrun_SCM_handler(
    m::ReferenceModel,
    tmpdir::String,
    u::Array{FT, 1},
    u_names::Array{String, 1},
    param_map::ParameterMap,
) where {FT<:AbstractFloat}Run a case using a set of parameters u_names with values u, and return directory pointing to where data is stored for simulation run.
Inputs:
- m :: Reference model
- tmpdir :: Temporary directory to store simulation results in
- u :: Values of parameters to be used in simulations.
- u_names       :: SCM names for parameters u.
- param_map:: A mapping operator to define relations between parameters. See- ParameterMapfor details.
Outputs:
- output_dir :: directory containing output data from the SCM run.
- model_error :: Boolean specifying whether the simulation failed.
run_SCM_handler(
    case_name::String,
    out_dir::String;
    u::Vector{FT},
    u_names::Vector{String},
    param_map::ParameterMap,
    namelist::Dict,
    uuid::String = "01",
    les::Union{NamedTuple, String} = nothing,
)Run a TurbulenceConvection.jl case and return directory pointing to where data is stored for simulation run.
Inputs:
- case_name :: case name
- out_dir :: Directory to store simulation results in.
Optional Inputs:
- u :: Values of parameters to be used in simulations.
- u_names       :: SCM names for parameters u.
- param_map:: A mapping operator to define relations between parameters. See- ParameterMapfor details.
- namelist :: namelist to use for simulation.
- uuid :: uuid of SCM run
- les :: path to LES stats file, or NamedTuple with keywords {forcingmodel, month, experiment, cfsitenumber} needed to specify path.
Outputs:
- output_dir :: directory containing output data from the SCM run.
- model_error :: Boolean specifying whether the simulation failed.
CalibrateEDMF.TurbulenceConvectionUtils.create_parameter_vectors — Functioncreate_parameter_vectors(u_names, u, param_map, namelist)Given vector of parameter names and corresponding values, combine any vector components into single parameter vectors for input into SCM.
Inputs:
- u_names:: SCM names for parameters- u, which may contain vector components.
- u:: Values of parameters to be used in simulations, which may contain vector components.
- param_map:: A mapping to a reduced parameter set. See- ParameterMapfor details.
- namelist:: The parameter namelist for TurbulenceConvection.jl
Outputs:
- u_names_out:: SCM names for parameters- u.
- u_out:: Values of parameters to be used in simulations.
CalibrateEDMF.TurbulenceConvectionUtils.generate_scm_input — Functiongenerate_scm_input(model_evaluators, iteration, outdir_path = pwd(), batch_indices = nothing)Writes to file a set of ModelEvaluator used to initialize SCM evaluations at different parameter vectors, as well as their assigned version, which is on the form i<iteration>_e<ensemble_index>.
CalibrateEDMF.TurbulenceConvectionUtils.parse_version_inds — Functionparse_version_inds(version)Given version = "ix_ey", return the iteration index x and ensemble index y as integers.
CalibrateEDMF.TurbulenceConvectionUtils.precondition — Functionprecondition(
    param::Vector{FT},
    priors,
    param_map::ParameterMap,
    ref_models::Vector{ReferenceModel},
    ref_stats::ReferenceStatistics;
    counter::Integer = 0,
    max_counter::Integer = 10,
) where {FT <: Real}Substitute parameter vector param by a parameter vector drawn from the same prior, conditioned on the forward model being stable.
Inputs:
- param:: A parameter vector that may possibly result in unstable forward model evaluations (in unconstrained space).
- priors:: Priors from which the parameters were drawn.
- param_map:: A mapping to a reduced parameter set. See- ParameterMapfor details.
- ref_models:: Vector of ReferenceModels to check stability for.
- ref_stats:: ReferenceStatistics of the ReferenceModels.
- counter:: Accumulator tracking number of recursive calls to preconditioner.
- max_counter:: Maximum number of recursive calls to the preconditioner.
Outputs:
- new_param:: A new parameter vector drawn from the prior, conditioned on simulations being stable (in unconstrained space).
precondition(ME::ModelEvaluator, priors)Substitute the parameter vector of a ModelEvaluator by another one drawn from the given priors, conditioned on the forward model being stable.
Inputs:
- ME :: A ModelEvaluator.
- priors :: Priors from which the parameters were drawn.
Outputs:
- A preconditioned ModelEvaluator.
CalibrateEDMF.TurbulenceConvectionUtils.get_gcm_les_uuid — Functionget_gcm_les_uuid(cfsite_number; [forcing_model::String, month, experiment])Generate unique and self-describing uuid given information about a GCM-driven LES simulation from [1].
Examples
julia> get_gcm_les_uuid(1; forcing_model = "HadGEM2-A", month = 7, experiment = "amip")
"1_HadGEM2-A_07_amip"CalibrateEDMF.TurbulenceConvectionUtils.save_tc_data — Functionsave_tc_data(cases, outdir_path, sim_dirs, version, suffix)Save full TC.jl output in <results_folder>/timeseries.<suffix>/iter_<iteration>/Output.<case>.<case_id>_<ens_i>.
Behavior of this function is specified in the output config, i.e. config["output"]. If the flag save_tc_output is set to true, TC.jl output is saved, and if additionally a list of iterations is specified in save_tc_iterations, only these EKP iterations are saved.
Arguments:
- config:: The calibration config dictionary. To save TC.jl output, set- save_tc_outputto- truein the output config. To only save specific EKP iterations, specify these in a vector in- save_tc_iterationsin the output config.
- outdir_path:: The results- results_folderpath
- sim_dirs:: List of (temporary) directories where raw simulation output is initially saved to
- version:: An identifier for the current iteration and ensemble index
- suffix:: Case set identifier; is either "train" or "validation".