5.13.3-dev0
mHM
The mesoscale Hydrological Model
Loading...
Searching...
No Matches
mo_objective_function Module Reference

Objective Functions for Optimization of mHM. More...

Functions/Subroutines

real(dp) function, public objective (parameterset, eval, arg1, arg2, arg3)
 Wrapper for objective functions.
 
real(dp) function objective_sm_kge_catchment_avg (parameterset, eval)
 Wrapper for objective functions.
 
real(dp) function, dimension(6) objective_q_et_tws_kge_catchment_avg (parameterset, eval)
 Objective function for et, tws and q.
 
subroutine init_indexarray_for_opti_data (domainmeta, optidataoption, noptidomains, opti_domain_indices)
 creates an index array of the inidices of the domains eval should MPI process.
 
real(dp) function objective_sm_corr (parameterset, eval)
 Objective function for soil moisture.
 
real(dp) function objective_sm_pd (parameterset, eval)
 Objective function for soil moisture.
 
real(dp) function objective_sm_sse_standard_score (parameterset, eval)
 Objective function for soil moisture.
 
real(dp) function objective_kge_q_rmse_tws (parameterset, eval)
 Objective function of KGE for runoff and RMSE for domain_avg TWS (standarized scores)
 
real(dp) function objective_neutrons_kge_catchment_avg (parameterset, eval)
 Objective function for neutrons.
 
real(dp) function objective_et_kge_catchment_avg (parameterset, eval)
 Objective function for evpotranspirstion (et).
 
real(dp) function objective_kge_q_sm_corr (parameterset, eval)
 Objective function of KGE for runoff and correlation for SM.
 
real(dp) function objective_kge_q_et (parameterset, eval)
 Objective function of KGE for runoff and KGE for ET.
 
real(dp) function objective_kge_q_bfi (parameterset, eval)
 Objective function of KGE for runoff and BFI absulute difference.
 
real(dp) function objective_kge_q_rmse_et (parameterset, eval)
 Objective function of KGE for runoff and RMSE for domain_avg ET (standarized scores)
 
subroutine create_domain_avg_tws (idomain, twsoptisim, tws_catch_avg_domain, tws_opti_catch_avg_domain, mask_times)
 
subroutine create_domain_avg_et (idomain, etoptisim, et_catch_avg_domain, et_opti_catch_avg_domain, mask_times)
 
subroutine convert_tws_to_twsa (twsoptisim, l1_twsaobs, twsaoptisim)
 

Detailed Description

Objective Functions for Optimization of mHM.

This module provides a wrapper for several objective functions used to optimize mHM against various variables. If the objective is only regarding runoff move it to mRM/mo_mrm_objective_function_runoff.f90. If it contains besides runoff another variable like TWS implement it here.

All the objective functions are supposed to be minimized!

  • (10) SO: SM: 1.0 - KGE of catchment average soilmoisture
  • (11) SO: SM: 1.0 - Pattern dissimilarity (PD) of spatially distributed soil moisture
  • (12) SO: SM: Sum of squared errors (SSE) of spatially distributed standard score (normalization) of soil moisture
  • (13) SO: SM: 1.0 - average temporal correlation of spatially distributed soil moisture
  • (15) SO: Q + TWS: [1.0-KGE(Q)]*RMSE(domain_avg_TWS) - objective function using Q and domain average (standard score) TWS
  • (17) SO: N: 1.0 - KGE of spatio-temporal neutron data, catchment-average
  • (27) SO: ET: 1.0 - KGE of catchment average evapotranspiration
    Changelog

Oldrich Rakovec Oct 2015

  • added obj. func. 15 (objective_kge_q_rmse_tws) and extract_domain_avg_tws routine, former basin_avg

Robert Schweppe Jun 2018

  • refactoring and reformatting
    Authors
    Juliane Mai
    Date
    Dec 2012

Function/Subroutine Documentation

◆ convert_tws_to_twsa()

subroutine mo_objective_function::convert_tws_to_twsa ( type(optidata_sim), intent(in) twsoptisim,
type(optidata), intent(in) l1_twsaobs,
type(optidata_sim), intent(inout) twsaoptisim )
private

Definition at line 2713 of file mo_objective_function.F90.

Referenced by objective_q_et_tws_kge_catchment_avg().

Here is the caller graph for this function:

◆ create_domain_avg_et()

subroutine mo_objective_function::create_domain_avg_et ( integer(i4), intent(in) idomain,
type(optidata_sim), dimension(:), intent(in) etoptisim,
real(dp), dimension(:), intent(out), allocatable et_catch_avg_domain,
real(dp), dimension(:), intent(out), allocatable et_opti_catch_avg_domain,
logical, dimension(:), intent(out), allocatable mask_times )
private

Definition at line 2658 of file mo_objective_function.F90.

References mo_global_variables::l1_etobs, and mo_common_constants::nodata_dp.

Referenced by objective_et_kge_catchment_avg().

Here is the caller graph for this function:

◆ create_domain_avg_tws()

subroutine mo_objective_function::create_domain_avg_tws ( integer(i4), intent(in) idomain,
type(optidata_sim), dimension(:), intent(in) twsoptisim,
real(dp), dimension(:), intent(out), allocatable tws_catch_avg_domain,
real(dp), dimension(:), intent(out), allocatable tws_opti_catch_avg_domain,
logical, dimension(:), intent(out), allocatable mask_times )
private

Definition at line 2603 of file mo_objective_function.F90.

References mo_global_variables::l1_twsaobs, and mo_common_constants::nodata_dp.

Referenced by objective_kge_q_rmse_tws().

Here is the caller graph for this function:

◆ init_indexarray_for_opti_data()

subroutine mo_objective_function::init_indexarray_for_opti_data ( type(domain_meta), intent(in) domainmeta,
integer(i4), intent(in) optidataoption,
integer(i4), intent(out) noptidomains,
integer(i4), dimension(:), intent(out), allocatable opti_domain_indices )
private

creates an index array of the inidices of the domains eval should MPI process.

The data type domainMeta contains an array optidata of size domainMetanDomains, telling us, which domains should be optimized with which opti_data. This subroutine splits all domains assigned to a process and returns an index list corresponding to the value of domainMetaoptidata.

The index array opti_domain_indices can then be passed as an optional argument to the eval subroutine. The eval then instead of using loops over all domains only uses the passed indices.

This subroutine also returns the size of that array since it helps with the calculations of the optimization in the end.

Authors
Maren Kaluza
Date
July 2019
Parameters
[in]domainmetameta data for all domains assigned to that process
[in]optidataoptionwhich opti data should be used in the eval called after calling this subroutine
[out]noptidomainsnumber of domains that will be optimized in the following eval call
[out]opti_domain_indicesthe indices of the domains that are to be processed in the following eval call

domain loop counter

Definition at line 900 of file mo_objective_function.F90.

Referenced by objective_q_et_tws_kge_catchment_avg().

Here is the caller graph for this function:

◆ objective()

real(dp) function, public mo_objective_function::objective ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval,
real(dp), intent(in), optional arg1,
real(dp), intent(out), optional arg2,
real(dp), intent(out), optional arg3 )

Wrapper for objective functions.

The functions selects the objective function case defined in a namelist, i.e. the global variable opti_function. It return the objective function value for a specific parameter set.

Parameters
[in]REAL(dp), DIMENSION(:) :: parameterset
[in]procedure(eval_interface) :: eval
[in]real(dp), optional :: arg1
[out]real(dp), optional :: arg2
[out]real(dp), optional :: arg3
Returns
real(dp) :: objective — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Juliane Mai
Date
Dec 2012

Definition at line 91 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_common_constants::nodata_dp, objective(), objective_et_kge_catchment_avg(), objective_kge_q_bfi(), objective_kge_q_et(), objective_kge_q_rmse_et(), objective_kge_q_rmse_tws(), objective_kge_q_sm_corr(), objective_neutrons_kge_catchment_avg(), objective_q_et_tws_kge_catchment_avg(), objective_sm_corr(), objective_sm_kge_catchment_avg(), objective_sm_pd(), objective_sm_sse_standard_score(), and mo_common_mhm_mrm_variables::opti_function.

Referenced by mo_mhm_interface::mhm_interface_run_optimization(), and objective().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_et_kge_catchment_avg()

real(dp) function mo_objective_function::objective_et_kge_catchment_avg ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function for evpotranspirstion (et).

The objective function only depends on a parameter vector. The model will be called with that parameter vector and the model output is subsequently compared to observed data. Therefore, the Kling-Gupta model efficiency \( KGE \) of the catchment average evapotranspiration (et) is calculated

\[ KGE = 1.0 - \sqrt{( (1-r)^2 + (1-\alpha)^2 + (1-\beta)^2 )} \]

where \( r \) = Pearson product-moment correlation coefficient \( \alpha \) = ratio of simulated mean to observed mean SM \( \beta \) = ratio of similated standard deviation to observed standard deviation is calculated and the objective function for a given domain \( i \) is

\[ \phi_{i} = 1.0 - KGE_{i} \]

\( \phi_{i} \) is the objective since we always apply minimization methods. The minimal value of \( \phi_{i} \) is 0 for the optimal KGE of 1.0. Finally, the overall objective function value \( OF \) is estimated based on the power-6 norm to combine the \( \phi_{i} \) from all domains \( N \).

\[ OF = \sqrt[6]{\sum((1.0 - KGE_{i})/N)^6 }. \]

The observed data L1_et, L1_et_mask are global in this module.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_et_kge_catchment_avg — objective function value (which will be e.g. minimized by an optimization routine)
Authors
Johannes Brenner
Date
Feb 2017

simulated et

Definition at line 1765 of file mo_objective_function.F90.

References create_domain_avg_et(), mo_common_variables::domainmeta, mo_common_constants::nodata_dp, and objective_et_kge_catchment_avg().

Referenced by objective(), and objective_et_kge_catchment_avg().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_kge_q_bfi()

real(dp) function mo_objective_function::objective_kge_q_bfi ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function of KGE for runoff and BFI absulute difference.

Objective function of KGE for runoff and KGE for ET. Further details can be found in the documentation of objective functions '14 - objective_multiple_gauges_kge_power6'.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_kge_q_BFI — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Sebastian Müller
Date
Apr 2022

baseflow index for each domain

Definition at line 2221 of file mo_objective_function.F90.

References mo_global_variables::bfi_obs, mo_common_variables::domainmeta, mo_mrm_objective_function_runoff::extract_runoff(), mo_common_variables::level1, and objective_kge_q_bfi().

Referenced by objective(), and objective_kge_q_bfi().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_kge_q_et()

real(dp) function mo_objective_function::objective_kge_q_et ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function of KGE for runoff and KGE for ET.

Objective function of KGE for runoff and KGE for ET. Further details can be found in the documentation of objective functions '14 - objective_multiple_gauges_kge_power6'.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_kge_q_et — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Johannes Brenner
Date
July 2017

simulated et

Definition at line 2045 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_mrm_objective_function_runoff::extract_runoff(), mo_global_variables::l1_etobs, mo_common_variables::level1, and objective_kge_q_et().

Referenced by objective(), and objective_kge_q_et().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_kge_q_rmse_et()

real(dp) function mo_objective_function::objective_kge_q_rmse_et ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function of KGE for runoff and RMSE for domain_avg ET (standarized scores)

Objective function of KGE for runoff and RMSE for domain_avg ET (standarized scores)

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_kge_q_rmse_et — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Johannes Brenner
Date
July 2017

simulated et

Definition at line 2353 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_common_constants::eps_dp, mo_common_mhm_mrm_variables::evalper, mo_mrm_objective_function_runoff::extract_runoff(), mo_global_variables::l1_etobs, mo_common_constants::nodata_dp, and objective_kge_q_rmse_et().

Referenced by objective(), and objective_kge_q_rmse_et().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_kge_q_rmse_tws()

real(dp) function mo_objective_function::objective_kge_q_rmse_tws ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function of KGE for runoff and RMSE for domain_avg TWS (standarized scores)

Objective function of KGE for runoff and RMSE for domain_avg TWS (standarized scores)

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_kge_q_rmse_tws — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Oldrich Rakovec, Rohini Kumar
Date
Oct. 2015

simulated tws

Definition at line 1383 of file mo_objective_function.F90.

References create_domain_avg_tws(), mo_common_variables::domainmeta, mo_common_constants::eps_dp, mo_common_mhm_mrm_variables::evalper, mo_mrm_objective_function_runoff::extract_runoff(), mo_global_variables::l1_twsaobs, mo_common_constants::nodata_dp, and objective_kge_q_rmse_tws().

Referenced by objective(), and objective_kge_q_rmse_tws().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_kge_q_sm_corr()

real(dp) function mo_objective_function::objective_kge_q_sm_corr ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function of KGE for runoff and correlation for SM.

Objective function of KGE for runoff and SSE for soil moisture (standarized scores). Further details can be found in the documentation of objective functions '14 - objective_multiple_gauges_kge_power6' and '13 - objective_sm_corr'.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_kge_q_sse_sm — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Matthias Zink
Date
Mar. 2017

Definition at line 1866 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_mrm_objective_function_runoff::extract_runoff(), mo_global_variables::l1_smobs, mo_common_variables::level1, mo_mrm_objective_function_runoff::objective_kge(), and objective_kge_q_sm_corr().

Referenced by objective(), and objective_kge_q_sm_corr().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_neutrons_kge_catchment_avg()

real(dp) function mo_objective_function::objective_neutrons_kge_catchment_avg ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function for neutrons.

The objective function only depends on a parameter vector. The model will be called with that parameter vector and the model output is subsequently compared to observed data. Therefore, the Kling-Gupta model efficiency \( KGE \) of the catchment average neutrons (N) is calculated

\[ KGE = 1.0 - \sqrt{( (1-r)^2 + (1-\alpha)^2 + (1-\beta)^2 )} \]

where \( r \) = Pearson product-moment CORRELATION coefficient \( \alpha \) = ratio of simulated mean to observed mean SM \( \beta \) = ratio of similated standard deviation to observed standard deviation is calculated and the objective function for a given domain \( i \) is

\[ \phi_{i} = 1.0 - KGE_{i} \]

\( \phi_{i} \) is the objective since we always apply minimization methods. The minimal value of \( \phi_{i} \) is 0 for the optimal KGE of 1.0. Finally, the overall objective function value \( OF \) is estimated based on the power-6 norm to combine the \( \phi_{i} \) from all domains \( N \).

\[ OF = \sqrt[6]{\sum((1.0 - KGE_{i})/N)^6 }. \]

The observed data L1_neutronsdata, L1_neutronsdata_mask are global in this module.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_neutrons_kge_catchment_avg — objective function value (which will be e.g. minimized by an optimization routine)
Authors
Martin Schroen
Date
Jun 2015

Definition at line 1619 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_global_variables::l1_neutronsobs, mo_common_constants::nodata_dp, and objective_neutrons_kge_catchment_avg().

Referenced by objective(), and objective_neutrons_kge_catchment_avg().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_q_et_tws_kge_catchment_avg()

real(dp) function, dimension(6) mo_objective_function::objective_q_et_tws_kge_catchment_avg ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function for et, tws and q.

The feature of this objective function is the separation of the eval call into four calls, each with another index list. The subroutine eval then only uses the indices from that index list internally instead of having loops over all domains. The integer array domainMetaoptidata decides which indices to use. Therefore the array is split into disjunct subsets, and, if chosen wisely in the namelist, also covers all domains.

With this the eval calls sum up in a way that for each domain eval was called at most once, but for different opti_data.

Authors
Maren Kaluza
Date
July 2019
Parameters
[in]parametersetthe parameterset passed to the eval subroutine
[in]evalthe eval subroutine called by this objective function
Returns
the return value of the objective function. In this case it is an array to provide the possibility to weight the outcome accordingly

domain loop counter

counter for short loops

for sixth root

modelled runoff for a given parameter set

number of all gauges, aquired via runoff

aggregated simulated runoff

measured runoff

mask for measured runoff

kge_q(nGaugesTotal)

gauges counter

number of q domains

number of et domains

number of tws domains

number of TWS and ET domains (providing both)

index array of ET domains

index array of TWS domains

index array of TWS and ET domains (providing both)

index array of ET domains

simulated et

simulated tws

simulated twsa (anomaly)

Definition at line 642 of file mo_objective_function.F90.

References convert_tws_to_twsa(), mo_common_variables::domainmeta, mo_mrm_objective_function_runoff::extract_runoff(), init_indexarray_for_opti_data(), mo_global_variables::l1_etobs, mo_global_variables::l1_twsaobs, mo_common_constants::nodata_dp, and objective_q_et_tws_kge_catchment_avg().

Referenced by objective(), and objective_q_et_tws_kge_catchment_avg().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_sm_corr()

real(dp) function mo_objective_function::objective_sm_corr ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function for soil moisture.

The objective function only depends on a parameter vector. The model will be called with that parameter vector and the model output is subsequently compared to observed data. Therefore the Pearson correlation between observed and modeled soil moisture on each grid cell \( j \) is compared

\[ r_j = r^2(SM_{obs}^j, SM_{sim}^j) \]

where \( r^2\) = Pearson correlation coefficient, \( SM_{obs} \) = observed soil moisture, \( SM_{sim} \) = simulated soil moisture. The observed data \( SM_{obs} \) are global in this module. The the correlation is spatially averaged as

\[ \phi_{i} = \frac{1}{K} \cdot \sum_{j=1}^K r_j \]

where \( K \) denotes the number of valid cells in the study domain. Finally, the overall objective function value \( OF \) is estimated based on the power-6 norm to combine the \( \phi_{i} \) from all domains \( N \).

\[ OF = \sqrt[6]{\sum((1.0 - \phi_{i})/N)^6 }. \]

The observed data L1_sm, L1_sm_mask are global in this module.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_sm_corr — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Matthias Zink
Date
March 2015

Definition at line 977 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_global_variables::l1_smobs, mo_common_variables::level1, and objective_sm_corr().

Referenced by objective(), and objective_sm_corr().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_sm_kge_catchment_avg()

real(dp) function mo_objective_function::objective_sm_kge_catchment_avg ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Wrapper for objective functions.

The functions sends the parameterset to the subprocess, receives the partial objective and calculates the final objective

Parameters
[in]REAL(dp), DIMENSION(:) :: parameterset
[in]procedure(eval_interface) :: eval
[in]real(dp), optional :: arg1
[out]real(dp), optional :: arg2
[out]real(dp), optional :: arg3
Returns
real(dp) :: objective — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Juliane Mai
Date
Dec 2012

Objective function for soil moisture.

The objective function only depends on a parameter vector. The model will be called with that parameter vector and the model output is subsequently compared to observed data. Therefore, the Kling-Gupta model efficiency \( KGE \) of the catchment average soil mloisture (SM) is calculated

\[ KGE = 1.0 - \sqrt{( (1-r)^2 + (1-\alpha)^2 + (1-\beta)^2 )} \]

where \( r \) = Pearson product-moment correlation coefficient \( \alpha \) = ratio of simulated mean to observed mean SM \( \beta \) = ratio of similated standard deviation to observed standard deviation is calculated and the objective function for a given domain \( i \) is

\[ \phi_{i} = 1.0 - KGE_{i} \]

\( \phi_{i} \) is the objective since we always apply minimization methods. The minimal value of \( \phi_{i} \) is 0 for the optimal KGE of 1.0. Finally, the overall objective function value \( OF \) is estimated based on the power-6 norm to combine the \( \phi_{i} \) from all domains \( N \).

\[ OF = \sqrt[6]{\sum((1.0 - KGE_{i})/N)^6 }. \]

The observed data L1_sm, L1_sm_mask are global in this module.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_sm_kge_catchment_avg — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Matthias Zink
Date
May 2015

Definition at line 498 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_global_variables::l1_smobs, mo_common_variables::level1, mo_common_constants::nodata_dp, and objective_sm_kge_catchment_avg().

Referenced by objective(), and objective_sm_kge_catchment_avg().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_sm_pd()

real(dp) function mo_objective_function::objective_sm_pd ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function for soil moisture.

The objective function only depends on a parameter vector. The model will be called with that parameter vector and the model output is subsequently compared to observed data. Therefore the Pattern Dissimilarity (PD) of observed and modeled soil moisture fields is calculated - aim: matching spatial patters

\[ E(t) = PD\left( SM_{obs}(t), SM_{sim}(t) \right) \]

where \( PD \) = pattern dissimilarity function, \( SM_{obs} \) = observed soil moisture, \( SM_{sim} \) = simulated soil moisture. \( E(t) \) = pattern dissimilarity at timestep \( t \). The the pattern dissimilaity (E) is spatially averaged as

\[ \phi_{i} = \frac{1}{T} \cdot \sum_{t=1}^T E_t \]

where \( T \) denotes the number of time steps. Finally, the overall objective function value \( OF \) is estimated based on the power-6 norm to combine the \( \phi_{i} \) from all domains \( N \).

\[ OF = \sqrt[6]{\sum((1.0 - \phi_{i})/N)^6 } . \]

The observed data L1_sm, L1_sm_mask are global in this module. The observed data L1_sm, L1_sm_mask are global in this module.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objecive_sm_pd — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Matthias Zink
Date
May 2015

Definition at line 1110 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_global_variables::l1_smobs, mo_common_variables::level1, mo_common_constants::nodata_dp, and objective_sm_pd().

Referenced by objective(), and objective_sm_pd().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ objective_sm_sse_standard_score()

real(dp) function mo_objective_function::objective_sm_sse_standard_score ( real(dp), dimension(:), intent(in) parameterset,
procedure(eval_interface), intent(in), pointer eval )
private

Objective function for soil moisture.

The objective function only depends on a parameter vector. The model will be called with that parameter vector and the model output is subsequently compared to observed data. Therefore the sum of squared errors (SSE) of the standard score of observed and modeled soil moisture is calculated. The standard score or normalization (anomaly) make the objctive function bias insensitive and basically the dynamics of the soil moisture is tried to capture by this obejective function.

\[ phi_i = \sum_{j=1}^K \{ standard\_score( SM_{obs}(j) )- standard\_score(SM_{sim}(j)) \}^2 \]

where \( standard\_score \) = standard score function, \( SM_{obs} \) = observed soil moisture, \( SM_{sim} \) = simulated soil moisture. \( K \) = valid elements in study domain. Finally, the overall objective function value \( OF \) is estimated based on the power-6 norm to combine the \( \phi_{i} \) from all domains \( N \).

\[ OF = \sqrt[6]{\sum(\phi_{i}/N)^6 }. \]

The observed data L1_sm, L1_sm_mask are global in this module.

Parameters
[in]real(dp), dimension(:) :: parameterset
[in]procedure(eval_interface) :: eval
Returns
real(dp) :: objective_sm_sse_standard_score — objective function value (which will be e.g. minimized by an optimization routine like DDS)
Authors
Matthias Zink
Date
March 2015

Definition at line 1264 of file mo_objective_function.F90.

References mo_common_variables::domainmeta, mo_global_variables::l1_smobs, mo_common_variables::level1, and objective_sm_sse_standard_score().

Referenced by objective(), and objective_sm_sse_standard_score().

Here is the call graph for this function:
Here is the caller graph for this function: