38 USE mo_kind,
ONLY : i4, dp
39 use mo_optimization_utils,
only : eval_interface
40 use mo_message,
only : message, error_message
48 PUBLIC :: objective_master, objective_subprocess
91 FUNCTION objective(parameterset, eval, arg1, arg2, arg3)
98 REAL(dp),
DIMENSION(:),
INTENT(IN) :: parameterset
100 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
102 real(dp),
optional,
intent(in) :: arg1
104 real(dp),
optional,
intent(out) :: arg2
106 real(dp),
optional,
intent(out) :: arg3
110 real(dp),
dimension(6) :: multiple_objective
113 if (
present(arg1) .or.
present(arg2) .or.
present(arg3))
then
114 call error_message(
"Error mo_objective_function: Received unexpected argument, check optimization settings")
118 if (
present(arg2))
then
121 if (
present(arg3))
then
164 call error_message(
"Error objective: opti_function not implemented yet.")
205 FUNCTION objective_master(parameterset, eval, arg1, arg2, arg3)
211 use mo_string_utils,
only : num2str
216 REAL(dp),
DIMENSION(:),
INTENT(IN) :: parameterset
218 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
220 real(dp),
optional,
intent(in) :: arg1
222 real(dp),
optional,
intent(out) :: arg2
224 real(dp),
optional,
intent(out) :: arg3
226 REAL(dp) :: objective_master
228 REAL(dp) :: partial_objective
230 real(dp),
dimension(6) :: multiple_partial_objective
232 real(dp),
dimension(6) :: multiple_master_objective
235 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
237 integer(i4) :: iproc, nproc
239 integer(i4) :: ierror
241 type(mpi_status) :: status
244 if (
present(arg1) .or.
present(arg2) .or.
present(arg3))
then
245 call error_message(
"Error mo_objective_function: Received unexpected argument, check optimization settings")
249 if (
present(arg2))
then
252 if (
present(arg3))
then
256 call distribute_parameterset(parameterset)
258 case (10 : 13, 17, 27 : 29)
259 call mpi_comm_size(
domainmeta%comMaster, nproc, ierror)
260 objective_master = 0.0_dp
261 do iproc = 1, nproc - 1
262 call mpi_recv(partial_objective, 1, mpi_double_precision, iproc, 0,
domainmeta%comMaster, status, ierror)
263 objective_master = objective_master + partial_objective
265 objective_master = objective_master**onesixth
268 call error_message(
"case 15, objective_kge_q_rmse_tws not implemented in parallel yet")
272 call message(
"case 30, objective_kge_q_rmse_et not implemented in parallel yet")
274 call mpi_comm_size(
domainmeta%comMaster, nproc, ierror)
275 objective_master = 0.0_dp
276 multiple_master_objective(:) = 0.0_dp
277 do iproc = 1, nproc - 1
278 call mpi_recv(multiple_partial_objective, 6, mpi_double_precision, iproc, 0,
domainmeta%comMaster, status, ierror)
279 multiple_master_objective = multiple_master_objective + multiple_partial_objective
281 objective_master = objective_master + &
282 (multiple_master_objective(1)+multiple_master_objective(2)+multiple_master_objective(3))
283 objective_master = (objective_master/multiple_master_objective(4))**onesixth
286 call error_message(
"Error objective_master: opti_function not implemented yet.")
291 call message(
' objective_sm_kge_catchment_avg = ', num2str(objective_master,
'(F9.5)'))
293 call message(
' objective_sm_pd = ', num2str(objective_master,
'(F9.5)'))
295 call message(
' objective_sm_sse_standard_score = ', num2str(objective_master,
'(E12.5)'))
297 call message(
' objective_sm_corr = ', num2str(objective_master,
'(F9.5)'))
299 call message(
' objective_neutrons_kge_catchment_avg = ', num2str(objective_master,
'(F9.5)'))
301 call message(
' objective_et_kge_catchment_avg = ', num2str(objective_master,
'(F9.5)'))
303 call message(
' objective_kge_q_sm_corr = ', num2str(objective_master,
'(F9.5)'))
305 call message(
' objective_kge_q_et = ', num2str(objective_master,
'(F9.5)'))
307 call message(
' objective_q_et_tws_kge_catchment_avg = ', num2str(objective_master,
'(F9.5)'))
309 call error_message(
"Error objective_master: opti_function not implemented yet, this part of the code should never execute.")
312 END FUNCTION objective_master
352 subroutine objective_subprocess(eval, arg1, arg2, arg3)
362 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
364 real(dp),
optional,
intent(in) :: arg1
366 real(dp),
optional,
intent(out) :: arg2
368 real(dp),
optional,
intent(out) :: arg3
370 REAL(dp) :: partial_objective
372 real(dp),
dimension(6) :: multiple_partial_objective
374 REAL(dp),
DIMENSION(:),
allocatable :: parameterset
376 integer(i4) :: ierror
378 type(mpi_status) :: status
380 logical :: do_obj_loop
383 call mpi_recv(do_obj_loop, 1, mpi_logical, 0, 0,
domainmeta%comMaster, status, ierror)
385 if (.not. do_obj_loop)
exit
387 if (
present(arg1) .or.
present(arg2) .or.
present(arg3))
then
388 call error_message(
"Error mo_objective_function: Received unexpected argument, check optimization settings")
392 if (
present(arg2))
then
395 if (
present(arg3))
then
398 call get_parameterset(parameterset)
415 call error_message(
"Error objective_subprocess: case 15 not supported with MPI.")
431 call error_message(
"Error objective_subprocess: case 30 not supported with MPI.")
435 call error_message(
"Error objective_subprocess: opti_function not implemented yet.")
439 case (10 : 13, 17, 27 : 29)
440 call mpi_send(partial_objective,1, mpi_double_precision,0,0,
domainmeta%comMaster,ierror)
442 call mpi_send(multiple_partial_objective, 6, mpi_double_precision,0,0,
domainmeta%comMaster,ierror)
444 call error_message(
"Error objective_subprocess: this part should not be executed -> error in the code.")
447 deallocate(parameterset)
450 END subroutine objective_subprocess
500 use mo_optimization_types,
only : optidata_sim
503 use mo_errormeasures,
only : kge
505 use mo_moment,
only : average
506 use mo_string_utils,
only : num2str
510 real(dp),
dimension(:),
intent(in) :: parameterset
512 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
517 integer(i4) :: idomain
523 integer(i4) :: n_time_steps
526 integer(i4) :: ncells1
529 real(dp) :: invalid_times
532 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
536 real(dp),
dimension(:),
allocatable :: sm_catch_avg_domain
539 real(dp),
dimension(:),
allocatable :: sm_opti_catch_avg_domain
541 type(optidata_sim),
dimension(:),
allocatable :: smoptisim
544 logical,
dimension(:),
allocatable :: mask_times
548 call eval(parameterset, smoptisim = smoptisim)
557 ncells1 =
level1(idomain)%ncells
560 allocate(mask_times(
size(smoptisim(idomain)%dataSim, dim = 2)))
561 allocate(sm_catch_avg_domain(
size(smoptisim(idomain)%dataSim, dim = 2)))
562 allocate(sm_opti_catch_avg_domain(
size(smoptisim(idomain)%dataSim, dim = 2)))
569 invalid_times = 0.0_dp
571 n_time_steps =
size(smoptisim(idomain)%dataSim, dim = 2)
572 do itime = 1, n_time_steps
576 if (count(
l1_smobs(idomain)%maskObs(:, itime)) .LE. (0.10_dp * real(ncells1, dp)))
then
577 invalid_times = invalid_times + 1.0_dp
578 mask_times(itime) = .false.
581 sm_catch_avg_domain(itime) = average(
l1_smobs(idomain)%dataObs(:, itime), mask =
l1_smobs(idomain)%maskObs(:, itime))
582 sm_opti_catch_avg_domain(itime) = average(smoptisim(idomain)%dataSim(:, itime), mask =
l1_smobs(idomain)%maskObs(:, itime))
586 if (invalid_times .GT. 0.5_dp)
then
587 call message(.LT.
' WARNING: objective_sm: Detected invalid timesteps ( 10 valid data points).')
588 call message(
' Fraction of invalid timesteps: ', &
589 num2str(invalid_times / real(n_time_steps, dp),
'(F4.2)'))
596 ((1.0_dp - kge(sm_catch_avg_domain, sm_opti_catch_avg_domain, mask = mask_times)) / &
597 real(
domainmeta%overallnumberofdomains, dp))**6
600 deallocate(mask_times)
601 deallocate(sm_catch_avg_domain)
602 deallocate(sm_opti_catch_avg_domain)
603 call smoptisim(idomain)%destroy()
605 deallocate(smoptisim)
644 use mo_optimization_types,
only : optidata_sim
648 use mo_errormeasures,
only : kge
649 use mo_moment,
only : average
650 use mo_string_utils,
only : num2str
656 real(dp),
dimension(:),
intent(in) :: parameterset
658 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
664 integer(i4) :: idomain
670 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
675 real(dp),
allocatable,
dimension(:, :) :: runoff
677 integer(i4) :: ngaugestotal
680 real(dp),
dimension(:),
allocatable :: runoff_agg
683 real(dp),
dimension(:),
allocatable :: runoff_obs
686 logical,
dimension(:),
allocatable :: runoff_obs_mask
692 integer(i4) :: gg, icell
695 integer(i4) :: nqdomains
698 integer(i4) :: netdomains
701 integer(i4) :: ntwsdomains
704 integer(i4) :: nettwsdomains
707 integer(i4),
dimension(:),
allocatable :: opti_domain_indices_et
710 integer(i4),
dimension(:),
allocatable :: opti_domain_indices_tws
713 integer(i4),
dimension(:),
allocatable :: opti_domain_indices_et_tws
716 integer(i4),
dimension(:),
allocatable :: opti_domain_indices_q
719 type(optidata_sim),
dimension(:),
allocatable :: etoptisim
722 type(optidata_sim),
dimension(:),
allocatable :: twsoptisim
725 type(optidata_sim),
dimension(:),
allocatable :: twsaoptisim
731 integer(i4) :: numberofsummands
747 if (nettwsdomains > 0)
then
751 call eval(parameterset, opti_domain_indices = opti_domain_indices_et_tws, &
752 twsoptisim = twsoptisim, etoptisim = etoptisim)
754 do i = 1,
size(opti_domain_indices_et_tws)
755 idomain = opti_domain_indices_et_tws(i)
757 do icell = 1,
size(
l1_etobs(idomain)%maskObs(:, :), dim = 1)
759 (1.0_dp - kge(
l1_etobs(idomain)%dataObs(icell, :), etoptisim(idomain)%dataSim(icell, :),&
760 mask =
l1_etobs(idomain)%maskObs(icell, :)))**6
761 numberofsummands = numberofsummands + 1
763 do icell = 1,
size(
l1_twsaobs(idomain)%maskObs(:, :), dim = 1)
764 kge_tws = kge_tws + &
765 (1.0_dp - kge(
l1_twsaobs(idomain)%dataObs(icell, :), twsaoptisim(idomain)%dataSim(icell, :),&
766 mask =
l1_twsaobs(idomain)%maskObs(icell, :)))**6
767 numberofsummands = numberofsummands + 1
770 call etoptisim(idomain)%destroy()
771 call twsoptisim(idomain)%destroy()
772 call twsaoptisim(idomain)%destroy()
774 deallocate(etoptisim)
775 deallocate(twsoptisim)
776 deallocate(twsaoptisim)
787 if (ntwsdomains > 0)
then
790 call eval(parameterset, opti_domain_indices = opti_domain_indices_tws, twsoptisim = twsoptisim)
792 do i = 1,
size(opti_domain_indices_tws)
793 idomain = opti_domain_indices_tws(i)
795 do icell = 1,
size(
l1_twsaobs(idomain)%maskObs(:, :), dim = 1)
796 kge_tws = kge_tws + &
797 (1.0_dp - kge(
l1_twsaobs(idomain)%dataObs(icell, :), twsaoptisim(idomain)%dataSim(icell, :),&
798 mask =
l1_twsaobs(idomain)%maskObs(icell, :)))**6
799 numberofsummands = numberofsummands + 1
801 call twsoptisim(idomain)%destroy()
803 deallocate(twsoptisim)
815 if (netdomains > 0)
then
817 call eval(parameterset, opti_domain_indices = opti_domain_indices_et, etoptisim = etoptisim)
819 do i = 1,
size(opti_domain_indices_et)
820 idomain = opti_domain_indices_et(i)
821 do icell = 1,
size(
l1_etobs(idomain)%maskObs(:, :), dim = 1)
823 (1.0_dp - kge(
l1_etobs(idomain)%dataObs(icell, :), etoptisim(idomain)%dataSim(icell, :),&
824 mask =
l1_etobs(idomain)%maskObs(icell, :)))**6
825 numberofsummands = numberofsummands + 1
827 call etoptisim(idomain)%destroy()
829 deallocate(etoptisim)
844 if (nqdomains > 0)
then
845 call eval(parameterset, opti_domain_indices = opti_domain_indices_q, runoff = runoff)
846 ngaugestotal =
size(runoff, dim = 2)
847 do gg = 1, ngaugestotal
850 call extract_runoff(gg, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
853 (1.0_dp - kge(runoff_obs, runoff_agg, mask = runoff_obs_mask))**6
854 numberofsummands = numberofsummands + 1
855 deallocate (runoff_agg, runoff_obs, runoff_obs_mask)
867 call message(
' objective_q_et_tws_kge_catchment_avg = ', &
905 integer(i4),
intent(in) :: optidataOption
907 integer(i4),
intent(out) :: nOptiDomains
909 integer(i4),
dimension(:),
allocatable,
intent(out) :: opti_domain_indices
912 integer(i4) :: iDomain, i
914 if (
allocated(opti_domain_indices))
deallocate(opti_domain_indices)
917 do idomain = 1, domainmeta%nDomains
918 if (domainmeta%optidata(idomain) == optidataoption) noptidomains = noptidomains + 1
921 if (noptidomains > 0)
then
922 allocate(opti_domain_indices(noptidomains))
924 do idomain = 1, domainmeta%nDomains
925 if (domainmeta%optidata(idomain) == optidataoption)
then
927 opti_domain_indices(i) = idomain
979 use mo_optimization_types,
only : optidata_sim
982 use mo_moment,
only : correlation
983 use mo_string_utils,
only : num2str
987 real(dp),
dimension(:),
intent(in) :: parameterset
989 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
994 integer(i4) :: idomain
1000 integer(i4) :: ncells1
1003 real(dp) :: invalid_cells
1006 real(dp) :: objective_sm_corr_domain
1009 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
1012 type(optidata_sim),
dimension(:),
allocatable :: smoptisim
1016 call eval(parameterset, smoptisim = smoptisim)
1025 objective_sm_corr_domain = 0.0_dp
1027 ncells1 =
level1(idomain)%ncells
1029 invalid_cells = 0.0_dp
1032 do icell = 1,
size(
l1_smobs(idomain)%maskObs(:, :), dim = 1)
1035 if (count(
l1_smobs(idomain)%maskObs(icell, :)) .LE. 0.10_dp * real(
size(
l1_smobs(idomain)%dataObs(:, :), dim = 2), dp))
then
1036 invalid_cells = invalid_cells + 1.0_dp
1039 objective_sm_corr_domain = objective_sm_corr_domain + &
1040 correlation(
l1_smobs(idomain)%dataObs(icell, :), smoptisim(idomain)%dataSim(icell, :), &
1041 mask =
l1_smobs(idomain)%maskObs(icell, :))
1045 if (invalid_cells .GT. 0.5_dp)
then
1046 call message(.LT.
' WARNING: objective_sm: Detected invalid cells in study area ( 10 valid data points).')
1047 call message(
' Fraction of invalid cells: ', &
1048 num2str(invalid_cells / real(ncells1, dp),
'(F4.2)'))
1055 ((1.0_dp - objective_sm_corr_domain / real(ncells1, dp)) / real(
domainmeta%overallNumberOfDomains, dp))**6
1112 use mo_optimization_types,
only : optidata_sim
1116 use mo_spatialsimilarity,
only : pd
1117 use mo_string_utils,
only : num2str
1121 real(dp),
dimension(:),
intent(in) :: parameterset
1123 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
1129 integer(i4) :: idomain
1132 integer(i4) :: itime
1135 integer(i4) :: nrows1, ncols1
1139 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
1143 real(dp),
dimension(:, :),
allocatable :: mat1, mat2
1146 real(dp),
dimension(:),
allocatable :: pd_time_series
1149 type(optidata_sim),
dimension(:),
allocatable :: smoptisim
1152 logical,
dimension(:, :),
allocatable :: mask1
1155 logical,
dimension(:, :),
allocatable :: mask_sm
1158 logical,
dimension(:),
allocatable :: mask_times
1162 call eval(parameterset, smoptisim = smoptisim)
1171 mask1 =
level1(idomain)%mask
1172 ncols1 =
level1(idomain)%ncols
1173 nrows1 =
level1(idomain)%nrows
1176 allocate(mask_times(
size(smoptisim(idomain)%dataSim, dim = 2)))
1177 allocate(pd_time_series(
size(smoptisim(idomain)%dataSim, dim = 2)))
1178 allocate(mat1(nrows1, ncols1))
1179 allocate(mat2(nrows1, ncols1))
1180 allocate(mask_sm(nrows1, ncols1))
1183 mask_times = .false.
1184 pd_time_series = 0.0_dp
1187 do itime = 1,
size(smoptisim(idomain)%dataSim, dim = 2)
1189 mat2 = unpack(smoptisim(idomain)%dataSim(:, itime), mask1,
nodata_dp)
1190 mask_sm = unpack(
l1_smobs(idomain)%maskObs(:, itime), mask1, .false.)
1191 pd_time_series = pd(mat1, mat2, mask = mask_sm, valid = mask_times(itime))
1194 if (count(mask_times) > 0_i4)
then
1197 ((1.0_dp - sum(pd_time_series, mask = mask_times) / real(count(mask_times), dp)) / &
1198 real(
domainmeta%overallnumberofdomains, dp))**6
1200 call error_message(
'***ERROR: mo_objective_funtion: objective_sm_pd: No soil moisture observations available!')
1204 deallocate(mask_times)
1205 deallocate(pd_time_series)
1209 call smoptisim(idomain)%destroy()
1211 deallocate(smoptisim)
1216 call message(
' objective_sm_pd = ', num2str(
objective_sm_pd,
'(F9.5)'))
1266 use mo_optimization_types,
only : optidata_sim
1268 use mo_errormeasures,
only : sse
1270 use mo_standard_score,
only : standard_score
1271 use mo_string_utils,
only : num2str
1275 real(dp),
dimension(:),
intent(in) :: parameterset
1277 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
1282 integer(i4) :: idomain
1285 integer(i4) :: icell
1288 integer(i4) :: ncells1
1291 real(dp) :: invalid_cells
1294 real(dp) :: objective_sm_sse_standard_score_domain
1297 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
1301 type(optidata_sim),
dimension(:),
allocatable :: smoptisim
1304 call eval(parameterset, smoptisim = smoptisim)
1313 objective_sm_sse_standard_score_domain = 0.0_dp
1315 ncells1 =
level1(idomain)%nCells
1317 invalid_cells = 0.0_dp
1319 do icell = 1,
size(
l1_smobs(idomain)%maskObs(:, :), dim = 1)
1322 if (count(
l1_smobs(idomain)%maskObs(icell, :)) .LE. (0.10_dp * real(
size(
l1_smobs(idomain)%dataObs, dim = 2), dp)))
then
1323 invalid_cells = invalid_cells + 1.0_dp
1326 objective_sm_sse_standard_score_domain = objective_sm_sse_standard_score_domain + &
1327 sse(standard_score(
l1_smobs(idomain)%dataObs(icell, :), mask =
l1_smobs(idomain)%maskObs(icell, :)), &
1328 standard_score(smoptisim(idomain)%dataSim(icell, :), mask =
l1_smobs(idomain)%maskObs(icell, :)), &
1329 mask =
l1_smobs(idomain)%maskObs(icell, :))
1334 if (invalid_cells .GT. 0.5_dp)
then
1335 call message(.LT.
' WARNING: objective_sm: Detected invalid cells in study area ( 10 valid data points).')
1336 call message(
' Fraction of invalid cells: ', &
1337 num2str(invalid_cells / real(ncells1, dp),
'(F4.2)'))
1343 (objective_sm_sse_standard_score_domain / real(
domainmeta%overallNumberOfDomains, dp))**6
1385 use mo_optimization_types,
only : optidata_sim
1390 use mo_errormeasures,
only : rmse
1391 use mo_julian,
only : caldat
1392 use mo_moment,
only : mean
1393 use mo_standard_score,
only : classified_standard_score
1394 use mo_string_utils,
only : num2str
1395 use mo_temporal_aggregation,
only : day2mon_average
1396 use mo_errormeasures,
only : kge
1401 real(dp),
dimension(:),
intent(in) :: parameterset
1403 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
1409 real(dp),
allocatable,
dimension(:, :) :: runoff
1412 type(optidata_sim),
dimension(:),
allocatable :: twsoptisim
1415 integer(i4) :: domainid, idomain, pp, mmm
1417 integer(i4) :: year, month, day
1419 real(dp),
dimension(domainMeta%nDomains) :: inittime
1422 real(dp),
dimension(:),
allocatable :: tws_catch_avg_domain
1425 real(dp),
dimension(:),
allocatable :: tws_opti_catch_avg_domain
1428 logical,
dimension(:),
allocatable :: tws_obs_mask
1431 integer(i4) :: nmonths
1434 integer(i4),
dimension(:),
allocatable :: month_classes
1437 real(dp),
DIMENSION(:),
allocatable :: tws_sim_m_anom, tws_obs_m_anom
1440 real(dp),
dimension(:),
allocatable :: rmse_tws
1443 real(dp) :: rmse_tws_avg, kge_q_avg
1445 integer(i4) :: ngaugestotal
1448 real(dp),
dimension(:),
allocatable :: runoff_agg
1451 real(dp),
dimension(:),
allocatable :: runoff_obs
1454 logical,
dimension(:),
allocatable :: runoff_obs_mask
1457 real(dp),
dimension(:),
allocatable :: kge_q
1464 call eval(parameterset, runoff = runoff, twsoptisim = twsoptisim)
1475 if (.not. (
l1_twsaobs(idomain)%timeStepInput == -2))
then
1476 call message(
'objective_kge_q_rmse_tws: current implementation of this subroutine only allows monthly timesteps')
1483 call create_domain_avg_tws(idomain, twsoptisim, tws_catch_avg_domain, tws_opti_catch_avg_domain, tws_obs_mask)
1486 if (count(tws_obs_mask) .lt. 12 * 2)
then
1487 call message(
'objective_kge_q_rmse_tws: Length of TWS data of domain ', trim(adjustl(num2str(domainid))), &
1488 ' less than 2 years: this is not recommended')
1492 inittime(idomain) = real(
evalper(idomain)%julStart, dp)
1495 call caldat(int(inittime(idomain)), yy = year, mm = month, dd = day)
1497 nmonths =
size(tws_obs_mask)
1499 allocate (month_classes(nmonths))
1500 allocate (tws_obs_m_anom(nmonths))
1501 allocate (tws_sim_m_anom(nmonths))
1503 month_classes(:) = 0
1510 month_classes(pp) = mmm
1511 if (mmm .LT. 12)
then
1519 tws_obs_m_anom = classified_standard_score(tws_opti_catch_avg_domain, month_classes, mask = tws_obs_mask)
1520 tws_sim_m_anom = classified_standard_score(tws_catch_avg_domain, month_classes, mask = tws_obs_mask)
1521 rmse_tws(idomain) = rmse(tws_sim_m_anom, tws_obs_m_anom, mask = tws_obs_mask)
1523 deallocate (month_classes)
1524 deallocate (tws_sim_m_anom)
1525 deallocate (tws_obs_m_anom)
1526 deallocate (tws_catch_avg_domain, tws_opti_catch_avg_domain, tws_obs_mask)
1528 call twsoptisim(idomain)%destroy()
1531 rmse_tws_avg = sum(rmse_tws(:), abs(rmse_tws -
nodata_dp) .gt.
eps_dp) / &
1533 deallocate(rmse_tws)
1534 deallocate(twsoptisim)
1540 ngaugestotal =
size(runoff, dim = 2)
1541 allocate(kge_q(ngaugestotal))
1544 do gg = 1, ngaugestotal
1547 call extract_runoff(gg, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
1550 pp = count(runoff_agg .ge. 0.0_dp)
1551 if (pp .lt. 365 * 2)
then
1552 call message(
'objective_kge_q_rmse_tws: The simulation at gauge ', trim(adjustl(num2str(gg))), &
1553 ' is not long enough. Please provide at least 730 days of data.')
1556 kge_q(gg) = kge(runoff_obs, runoff_agg, mask = runoff_obs_mask)
1557 deallocate (runoff_agg, runoff_obs, runoff_obs_mask)
1621 use mo_optimization_types,
only : optidata_sim
1624 use mo_errormeasures,
only : kge
1626 use mo_moment,
only : average
1627 use mo_string_utils,
only : num2str
1631 real(dp),
dimension(:),
intent(in) :: parameterset
1633 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
1638 integer(i4) :: idomain
1641 integer(i4) :: itime
1645 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
1649 real(dp),
dimension(:),
allocatable :: neutrons_catch_avg_domain
1652 real(dp),
dimension(:),
allocatable :: neutrons_opti_catch_avg_domain
1655 type(optidata_sim),
dimension(:),
allocatable :: neutronsoptisim
1658 logical,
dimension(:),
allocatable :: mask_times
1661 allocate(neutronsoptisim(
domainmeta%nDomains))
1662 call eval(parameterset, neutronsoptisim = neutronsoptisim)
1671 allocate(mask_times(
size(neutronsoptisim(idomain)%dataSim, dim = 2)))
1672 allocate(neutrons_catch_avg_domain(
size(neutronsoptisim(idomain)%dataSim, dim = 2)))
1673 allocate(neutrons_opti_catch_avg_domain(
size(neutronsoptisim(idomain)%dataSim, dim = 2)))
1678 neutrons_opti_catch_avg_domain =
nodata_dp
1681 do itime = 1,
size(neutronsoptisim(idomain)%dataSim, dim = 2)
1685 call message(
'WARNING: neutrons data at time ', num2str(itime,
'(I10)'),
' is empty.')
1688 mask_times(itime) = .false.
1691 neutrons_catch_avg_domain(itime) = average(
l1_neutronsobs(idomain)%dataObs(:, itime), &
1693 neutrons_opti_catch_avg_domain(itime) = average(neutronsoptisim(idomain)%dataSim(:, itime), &
1700 ((1.0_dp - kge(neutrons_catch_avg_domain, neutrons_opti_catch_avg_domain, mask = mask_times)) / &
1701 real(
domainmeta%overallnumberofdomains, dp))**6
1704 deallocate(mask_times)
1705 deallocate(neutrons_catch_avg_domain)
1706 deallocate(neutrons_opti_catch_avg_domain)
1708 call neutronsoptisim(idomain)%destroy()
1710 deallocate(neutronsoptisim)
1767 use mo_optimization_types,
only : optidata_sim
1770 use mo_errormeasures,
only : kge
1771 use mo_moment,
only : average
1772 use mo_string_utils,
only : num2str
1776 real(dp),
dimension(:),
intent(in) :: parameterset
1778 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
1783 integer(i4) ::idomain
1787 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
1791 real(dp),
dimension(:),
allocatable :: et_catch_avg_domain
1794 real(dp),
dimension(:),
allocatable :: et_opti_catch_avg_domain
1797 type(optidata_sim),
dimension(:),
allocatable :: etoptisim
1800 logical,
dimension(:),
allocatable :: mask_times
1803 call eval(parameterset, etoptisim = etoptisim)
1814 et_opti_catch_avg_domain, mask_times)
1819 ((1.0_dp - kge(et_catch_avg_domain, et_opti_catch_avg_domain, mask = mask_times)) / &
1820 real(
domainmeta%overallnumberofdomains, dp))**6
1823 deallocate(mask_times)
1824 deallocate(et_catch_avg_domain)
1825 deallocate(et_opti_catch_avg_domain)
1826 call etoptisim(idomain)%destroy()
1828 deallocate(etoptisim)
1868 use mo_optimization_types,
only : optidata_sim
1871 use mo_moment,
only : correlation
1872 use mo_string_utils,
only : num2str
1873 use mo_errormeasures,
only : kge
1878 real(dp),
dimension(:),
intent(in) :: parameterset
1880 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
1884 real(dp) :: objective_sm
1889 real(dp) :: invalid_cells
1893 real(dp),
allocatable,
dimension(:, :) :: runoff
1896 integer(i4) :: idomain
1899 integer(i4) :: icell
1902 integer(i4) :: ncells1
1905 real(dp) :: objective_sm_domain
1908 type(optidata_sim),
dimension(:),
allocatable :: smoptisim
1910 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
1915 integer(i4) :: ngaugestotal
1918 real(dp),
dimension(:),
allocatable :: runoff_agg
1921 real(dp),
dimension(:),
allocatable :: runoff_obs
1924 logical,
dimension(:),
allocatable :: runoff_obs_mask
1928 call eval(parameterset, runoff = runoff, smoptisim = smoptisim)
1935 objective_sm = 0.0_dp
1941 objective_sm_domain = 0.0_dp
1943 ncells1 =
level1(idomain)%nCells
1946 invalid_cells = 0.0_dp
1947 do icell = 1,
size(
l1_smobs(idomain)%maskObs(:, :), dim = 1)
1950 if (count(
l1_smobs(idomain)%maskObs(icell, :)) .LE. (0.10_dp * real(
size(
l1_smobs(idomain)%dataObs, dim = 2), dp)))
then
1951 invalid_cells = invalid_cells + 1.0_dp
1956 objective_sm_domain = objective_sm_domain + &
1957 correlation(
l1_smobs(idomain)%dataObs(icell, :), smoptisim(idomain)%dataSim(icell, :), &
1958 mask =
l1_smobs(idomain)%maskObs(icell, :))
1962 if (invalid_cells .GT. 0.5_dp)
then
1963 call message(.LT.
' WARNING: objective_sm: Detected invalid cells in study area ( 10 valid data points).')
1964 call message(
' Fraction of invalid cells: ', &
1965 num2str(invalid_cells / real(ncells1, dp),
'(F4.2)'))
1970 objective_sm = objective_sm + &
1971 ((1.0_dp - objective_sm_domain / real(ncells1, dp)) / real(
domainmeta%overallNumberOfDomains, dp))**6
1972 call smoptisim(idomain)%destroy()
1974 deallocate(smoptisim)
1977 objective_sm = objective_sm**onesixth
1983 ngaugestotal =
size(runoff, dim = 2)
1985 do gg = 1, ngaugestotal
1988 call extract_runoff(gg, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
1992 ((1.0_dp - kge(runoff_obs, runoff_agg, mask = runoff_obs_mask)) / real(ngaugestotal, dp))**6
1996 deallocate(runoff_agg, runoff_obs, runoff_obs_mask)
2047 use mo_optimization_types,
only : optidata_sim
2049 use mo_errormeasures,
only : kge
2051 use mo_string_utils,
only : num2str
2056 real(dp),
dimension(:),
intent(in) :: parameterset
2058 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
2062 real(dp) :: objective_et
2064 real(dp) :: objective_q
2067 real(dp) :: invalid_cells
2071 real(dp),
allocatable,
dimension(:, :) :: runoff
2074 integer(i4) :: idomain
2077 integer(i4) :: icell
2080 integer(i4) :: ncells1
2083 real(dp) :: objective_et_domain
2086 type(optidata_sim),
dimension(:),
allocatable :: etoptisim
2088 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
2093 integer(i4) :: ngaugestotal
2096 real(dp),
dimension(:),
allocatable :: runoff_agg
2099 real(dp),
dimension(:),
allocatable :: runoff_obs
2102 logical,
dimension(:),
allocatable :: runoff_obs_mask
2106 call eval(parameterset, runoff = runoff, etoptisim = etoptisim)
2113 objective_et = 0.0_dp
2119 objective_et_domain = 0.0_dp
2121 ncells1 =
level1(idomain)%nCells
2124 invalid_cells = 0.0_dp
2125 do icell = 1,
size(
l1_etobs(idomain)%maskObs, dim = 1)
2128 if (count(
l1_etobs(idomain)%maskObs(icell, :)) .LE. &
2129 (0.10_dp * real(
size(
l1_etobs(idomain)%dataObs(:, :), dim = 2), dp)))
then
2130 invalid_cells = invalid_cells + 1.0_dp
2135 objective_et_domain = objective_et_domain + &
2136 kge(
l1_etobs(idomain)%dataObs(icell, :), etoptisim(idomain)%dataSim(icell, :), &
2137 mask =
l1_etobs(idomain)%maskObs(icell, :))
2141 if (invalid_cells .GT. 0.5_dp)
then
2142 call message(.LT.
' WARNING: objective_et: Detected invalid cells in study area ( 10 valid data points).')
2143 call message(
' Fraction of invalid cells: ', &
2144 num2str(invalid_cells / real(ncells1, dp),
'(F4.2)'))
2149 objective_et = objective_et + &
2150 ((1.0_dp - objective_et_domain / real(ncells1, dp)) / real(
domainmeta%overallNumberOfDomains, dp))**6
2151 call etoptisim(idomain)%destroy()
2153 deallocate(etoptisim)
2156 objective_et = objective_et**onesixth
2161 objective_q = 0.0_dp
2162 ngaugestotal =
size(runoff, dim = 2)
2164 do gg = 1, ngaugestotal
2167 call extract_runoff(gg, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
2170 objective_q = objective_q + &
2171 ((1.0_dp - kge(runoff_obs, runoff_agg, mask = runoff_obs_mask)) / real(ngaugestotal, dp))**6
2175 deallocate(runoff_agg, runoff_obs, runoff_obs_mask)
2178 objective_q = objective_q**onesixth
2223 use mo_optimization_types,
only : optidata_sim
2225 use mo_errormeasures,
only : kge
2227 use mo_string_utils,
only : num2str
2232 real(dp),
dimension(:),
intent(in) :: parameterset
2234 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
2237 real(dp) :: objective_bfi
2238 real(dp) :: objective_q
2240 real(dp),
parameter :: onesixth = 1.0_dp / 6.0_dp
2243 real(dp) :: invalid_cells
2247 real(dp),
allocatable,
dimension(:, :) :: runoff
2250 integer(i4) :: idomain
2253 real(dp),
dimension(:),
allocatable :: bfi
2256 integer(i4) :: gg, i
2258 integer(i4) :: ngaugestotal
2261 integer(i4),
dimension(:),
allocatable :: domain_ids, domain_ids_pack
2264 real(dp),
dimension(:),
allocatable :: runoff_agg
2267 real(dp),
dimension(:),
allocatable :: runoff_obs
2270 logical,
dimension(:),
allocatable :: runoff_obs_mask
2274 call eval(parameterset, runoff = runoff, bfi = bfi)
2281 objective_bfi = 0.0_dp
2283 if ( any(
bfi_obs < 0.0_dp) )
then
2285 allocate(domain_ids_pack(count(
bfi_obs < 0.0_dp)))
2286 domain_ids = [(i, i=1,
size(domain_ids))]
2287 domain_ids_pack = pack(domain_ids, mask=(
bfi_obs < 0.0_dp))
2288 call error_message( &
2289 "objective_kge_q_BFI: missing BFI values for domain ", &
2290 trim(adjustl(num2str(domain_ids_pack(1)))) &
2296 objective_bfi = objective_bfi + abs(bfi(idomain) -
bfi_obs(idomain)) /
domainmeta%nDomains
2302 objective_q = 0.0_dp
2303 ngaugestotal =
size(runoff, dim = 2)
2305 do gg = 1, ngaugestotal
2308 call extract_runoff(gg, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
2311 objective_q = objective_q + &
2312 ((1.0_dp - kge(runoff_obs, runoff_agg, mask = runoff_obs_mask)) / real(ngaugestotal, dp))**6
2316 deallocate(runoff_agg, runoff_obs, runoff_obs_mask)
2319 objective_q = objective_q**onesixth
2355 use mo_optimization_types,
only : optidata_sim
2359 use mo_errormeasures,
only : rmse
2361 use mo_julian,
only : caldat
2362 use mo_moment,
only : average, mean
2363 use mo_standard_score,
only : classified_standard_score
2364 use mo_string_utils,
only : num2str
2365 use mo_temporal_aggregation,
only : day2mon_average
2366 use mo_errormeasures,
only : kge
2371 real(dp),
dimension(:),
intent(in) :: parameterset
2373 procedure(eval_interface),
INTENT(IN),
POINTER :: eval
2379 real(dp),
allocatable,
dimension(:, :) :: runoff
2382 type(optidata_sim),
dimension(:),
allocatable :: etoptisim
2385 integer(i4) :: itime
2388 integer(i4) :: idomain, pp, mmm
2390 integer(i4) :: year, month, day
2392 real(dp),
dimension(domainMeta%nDomains) :: inittime
2395 integer(i4) :: nmonths
2398 integer(i4),
dimension(:),
allocatable :: month_classes
2401 real(dp),
dimension(:),
allocatable :: et_sim_m, et_obs_m
2404 real(dp),
dimension(:),
allocatable :: et_sim_m_anom, et_obs_m_anom
2406 logical,
dimension(:),
allocatable :: et_obs_m_mask
2409 real(dp),
dimension(:),
allocatable :: rmse_et
2412 real(dp) :: rmse_et_avg, kge_q_avg
2415 real(dp),
dimension(:),
allocatable :: et_catch_avg_domain
2418 real(dp),
dimension(:),
allocatable :: et_opti_catch_avg_domain
2421 logical,
dimension(:),
allocatable :: mask_times
2424 real(dp),
dimension(:),
allocatable :: kge_q
2429 integer(i4) :: ngaugestotal
2432 real(dp),
dimension(:),
allocatable :: runoff_agg
2435 real(dp),
dimension(:),
allocatable :: runoff_obs
2438 logical,
dimension(:),
allocatable :: runoff_obs_mask
2444 call eval(parameterset, runoff = runoff, etoptisim = etoptisim)
2457 allocate(mask_times(
size(etoptisim(idomain)%dataSim, dim = 2)))
2458 allocate(et_catch_avg_domain(
size(etoptisim(idomain)%dataSim, dim = 2)))
2459 allocate(et_opti_catch_avg_domain(
size(etoptisim(idomain)%dataSim, dim = 2)))
2467 do itime = 1,
size(etoptisim(idomain)%dataSim, dim = 2)
2469 if (all(.NOT.
l1_etobs(idomain)%maskObs(:, itime)))
then
2473 mask_times(itime) = .false.
2477 et_catch_avg_domain(itime) = average(
l1_etobs(idomain)%dataObs(:, itime), &
2478 mask =
l1_etobs(idomain)%maskObs(:, itime))
2480 et_opti_catch_avg_domain(itime) = average(etoptisim(idomain)%dataSim(:, itime), &
2481 mask =
l1_etobs(idomain)%maskObs(:, itime))
2485 inittime(idomain) = real(
evalper(idomain)%julStart, dp)
2488 call caldat(int(inittime(idomain)), yy = year, mm = month, dd = day)
2491 select case(
l1_etobs(idomain)%timeStepInput)
2496 call day2mon_average(et_opti_catch_avg_domain, year, month, day, et_sim_m, misval =
nodata_dp)
2498 call day2mon_average(et_catch_avg_domain, year, month, day, et_obs_m, misval =
nodata_dp)
2503 allocate(et_sim_m(
size(etoptisim(idomain)%dataSim, dim = 2)))
2504 et_sim_m = et_opti_catch_avg_domain
2506 allocate(et_obs_m(
size(etoptisim(idomain)%dataSim, dim = 2)))
2507 et_obs_m = et_catch_avg_domain
2511 call error_message(
'***ERROR: objective_kge_q_rmse_et: time step of evapotranspiration yearly.')
2514 et_sim_m(:) = et_sim_m(:) - mean(et_sim_m(:))
2516 et_obs_m(:) = et_obs_m(:) - mean(et_obs_m(:))
2518 nmonths =
size(et_obs_m)
2519 allocate (month_classes(nmonths))
2520 allocate (et_obs_m_mask(nmonths))
2521 allocate (et_obs_m_anom(nmonths))
2522 allocate (et_sim_m_anom(nmonths))
2524 month_classes(:) = 0
2525 et_obs_m_mask(:) = .true.
2532 month_classes(pp) = mmm
2533 if (mmm .LT. 12)
then
2544 et_obs_m_anom = classified_standard_score(et_obs_m, month_classes, mask = et_obs_m_mask)
2545 et_sim_m_anom = classified_standard_score(et_sim_m, month_classes, mask = et_obs_m_mask)
2546 rmse_et(idomain) = rmse(et_sim_m_anom, et_obs_m_anom, mask = et_obs_m_mask)
2548 deallocate (month_classes)
2549 deallocate (et_obs_m)
2550 deallocate (et_sim_m)
2551 deallocate (et_obs_m_mask)
2552 deallocate (et_sim_m_anom)
2553 deallocate (et_obs_m_anom)
2556 call etoptisim(idomain)%destroy()
2559 rmse_et_avg = sum(rmse_et(:), abs(rmse_et -
nodata_dp) .gt.
eps_dp) / &
2562 deallocate(etoptisim)
2568 ngaugestotal =
size(runoff, dim = 2)
2569 allocate(kge_q(ngaugestotal))
2572 do gg = 1, ngaugestotal
2575 call extract_runoff(gg, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
2585 kge_q(gg) = kge(runoff_obs, runoff_agg, mask = runoff_obs_mask)
2586 deallocate (runoff_agg, runoff_obs, runoff_obs_mask)
2604 tws_opti_catch_avg_domain, mask_times)
2605 use mo_optimization_types,
only : optidata_sim
2608 use mo_moment,
only : average
2610 integer(i4),
intent(in) :: iDomain
2613 type(optidata_sim),
dimension(:),
intent(in) :: twsOptiSim
2616 real(dp),
dimension(:),
allocatable,
intent(out) :: tws_catch_avg_domain
2619 real(dp),
dimension(:),
allocatable,
intent(out) :: tws_opti_catch_avg_domain
2622 logical,
dimension(:),
allocatable,
intent(out) :: mask_times
2626 integer(i4) :: iTime
2629 allocate(mask_times(
size(twsoptisim(idomain)%dataSim, dim = 2)))
2630 allocate(tws_catch_avg_domain(
size(twsoptisim(idomain)%dataSim, dim = 2)))
2631 allocate(tws_opti_catch_avg_domain(
size(twsoptisim(idomain)%dataSim, dim = 2)))
2639 do itime = 1,
size(twsoptisim(idomain)%dataSim, dim = 2)
2642 if (all(.NOT.
l1_twsaobs(idomain)%maskObs(:, itime)))
then
2646 mask_times(itime) = .false.
2650 tws_catch_avg_domain(itime) = average(
l1_twsaobs(idomain)%dataObs(:, itime), &
2651 mask =
l1_twsaobs(idomain)%maskObs(:, itime))
2652 tws_opti_catch_avg_domain(itime) = average(twsoptisim(idomain)%dataSim(:, itime), &
2653 mask =
l1_twsaobs(idomain)%maskObs(:, itime))
2659 et_opti_catch_avg_domain, mask_times)
2660 use mo_optimization_types,
only : optidata_sim
2663 use mo_moment,
only : average
2665 integer(i4),
intent(in) :: iDomain
2668 type(optidata_sim),
dimension(:),
intent(in) :: etOptiSim
2671 real(dp),
dimension(:),
allocatable,
intent(out) :: et_catch_avg_domain
2674 real(dp),
dimension(:),
allocatable,
intent(out) :: et_opti_catch_avg_domain
2677 logical,
dimension(:),
allocatable,
intent(out) :: mask_times
2681 integer(i4) :: iTime
2684 allocate(mask_times(
size(etoptisim(idomain)%dataSim, dim = 2)))
2685 allocate(et_catch_avg_domain(
size(etoptisim(idomain)%dataSim, dim = 2)))
2686 allocate(et_opti_catch_avg_domain(
size(etoptisim(idomain)%dataSim, dim = 2)))
2694 do itime = 1,
size(etoptisim(idomain)%dataSim, dim = 2)
2697 if (all(.NOT.
l1_etobs(idomain)%maskObs(:, itime)))
then
2701 mask_times(itime) = .false.
2705 et_catch_avg_domain(itime) = average(
l1_etobs(idomain)%dataObs(:, itime), &
2706 mask =
l1_etobs(idomain)%maskObs(:, itime))
2707 et_opti_catch_avg_domain(itime) = average(etoptisim(idomain)%dataSim(:, itime), &
2708 mask =
l1_etobs(idomain)%maskObs(:, itime))
2714 use mo_optimization_types,
only : optidata_sim, optidata
2715 use mo_moment,
only : average
2717 type(optidata_sim),
intent(in) :: twsOptiSim
2719 type(optidata),
intent(in) :: L1_twsaObs
2721 type(optidata_sim),
intent(inout) :: twsaOptiSim
2724 integer(i4) :: iCell
2725 real(dp) :: twsa_av_cell
2727 allocate(twsaoptisim%dataSim(
size(twsoptisim%dataSim(:, :), dim = 1),
size(twsoptisim%dataSim(:, :), dim = 2)))
2729 do icell = 1,
size(twsoptisim%dataSim(:, :), dim = 1)
2730 twsa_av_cell = average(twsoptisim%dataSim(icell, :), mask = l1_twsaobs%maskObs(icell, :))
2731 twsaoptisim%dataSim(icell, :) = twsoptisim%dataSim(icell, :) - twsa_av_cell
Provides constants commonly used by mHM, mRM and MPR.
real(dp), parameter, public eps_dp
epsilon(1.0) in double precision
real(dp), parameter, public nodata_dp
Provides structures needed by mHM, mRM and/or mpr.
integer(i4), public opti_function
type(period), dimension(:), allocatable, public evalper
Provides common types needed by mHM, mRM and/or mpr.
Provides structures needed by mHM, mRM and/or mpr.
type(domain_meta), public domainmeta
type(grid), dimension(:), allocatable, target, public level1
Main global variables for mHM.
type(optidata), dimension(:), allocatable, public l1_twsaobs
this stores L1_tws, the mask, the directory of the observerd data, and the timestepInput of the simul...
type(optidata), dimension(:), allocatable, public l1_neutronsobs
type(optidata), dimension(:), allocatable, public l1_smobs
type(optidata), dimension(:), allocatable, public l1_etobs
real(dp), dimension(:), allocatable, public bfi_obs
given base-flow index per domain
Objective Functions for Optimization of mHM/mRM against runoff.
subroutine, public extract_runoff(gaugeid, runoff, runoff_agg, runoff_obs, runoff_obs_mask)
extracts runoff data from global variables
real(dp) function objective_kge(parameterset, eval)
Objective function of KGE.
Objective Functions for Optimization of mHM.
real(dp) function objective_sm_corr(parameterset, eval)
Objective function for soil moisture.
real(dp) function, public objective(parameterset, eval, arg1, arg2, arg3)
Wrapper for objective functions.
subroutine create_domain_avg_et(idomain, etoptisim, et_catch_avg_domain, et_opti_catch_avg_domain, mask_times)
real(dp) function objective_neutrons_kge_catchment_avg(parameterset, eval)
Objective function for neutrons.
real(dp) function objective_sm_pd(parameterset, eval)
Objective function for soil moisture.
subroutine convert_tws_to_twsa(twsoptisim, l1_twsaobs, twsaoptisim)
real(dp) function objective_kge_q_rmse_et(parameterset, eval)
Objective function of KGE for runoff and RMSE for domain_avg ET (standarized scores)
real(dp) function objective_kge_q_sm_corr(parameterset, eval)
Objective function of KGE for runoff and correlation for SM.
subroutine create_domain_avg_tws(idomain, twsoptisim, tws_catch_avg_domain, tws_opti_catch_avg_domain, mask_times)
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_kge_q_bfi(parameterset, eval)
Objective function of KGE for runoff and BFI absulute difference.
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, dimension(6) objective_q_et_tws_kge_catchment_avg(parameterset, eval)
Objective function for et, tws and q.
real(dp) function objective_sm_kge_catchment_avg(parameterset, eval)
Wrapper for objective functions.
real(dp) function objective_et_kge_catchment_avg(parameterset, eval)
Objective function for evpotranspirstion (et).
real(dp) function objective_kge_q_et(parameterset, eval)
Objective function of KGE for runoff and KGE for ET.
DOMAIN general description.