21 vw_w,
tu_su,
uu_advection,
uu_viscosity,
wu_u,
tv_sv,
vv_advection,
vv_viscosity,
wv_v,
tw_sw,
ww_advection, &
49 integer :: current_index, total_number_published_fields
69 current_index=current_index+1
75 current_index=current_index+1
81 current_index=current_index+1
87 current_index=current_index+1
93 current_index=current_index+1
99 current_index=current_index+1
105 current_index=current_index+1
111 current_index=current_index+1
118 type(model_state_type),
target,
intent(inout) :: current_state
134 type(model_state_type),
target,
intent(inout) :: current_state
137 if (current_state%first_timestep_column)
then 147 if (.not. current_state%halo_column)
then 163 type(model_state_type),
target,
intent(inout) :: current_state
184 if (
allocated(
uw_w))
deallocate(
uw_w)
185 if (
allocated(
vw_w))
deallocate(
vw_w)
190 if (
allocated(
wu_u))
deallocate(
wu_u)
194 if (
allocated(
wv_v))
deallocate(
wv_v)
252 if (
allocated(
w_qt))
deallocate(
w_qt)
375 type(model_state_type),
target,
intent(inout) :: current_state
376 character(len=*),
intent(in) :: name
377 type(component_field_information_type),
intent(out) :: field_information
380 field_information%field_type=component_array_field_type
381 field_information%data_type=component_double_data_type
385 field_information%number_dimensions=1
386 field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
403 field_information%number_dimensions=2
404 field_information%dimension_sizes(1)=current_state%local_grid%size(z_index)
405 field_information%dimension_sizes(2)=current_state%number_q_fields
413 type(model_state_type),
target,
intent(inout) :: current_state
415 wmfcrit=options_get_real(current_state%options_database,
"wmfcrit")
422 type(model_state_type),
target,
intent(inout) :: current_state
424 integer :: column_size
425 logical :: us_qt_enabled, u_qt_advection_enabled, u_qt_viscosity_diffusion_enabled, &
426 wu_qt_enabled, vs_qt_enabled, v_qt_advection_enabled, &
427 v_qt_viscosity_diffusion_enabled, wv_qt_enabled, w_qt_enabled, ws_qt_enabled, &
428 w_qt_advection_enabled, w_qt_viscosity_diffusion_enabled, w_qt_buoyancy_enabled, ww_qt_enabled, &
429 qt_qt_enabled, sqt_qt_enabled, qt_qt_advection_enabled, &
430 qt_qt_diffusion_enabled, wqt_qt_enabled
432 column_size=current_state%local_grid%size(z_index)
434 us_qt_enabled=current_state%u%active .and. current_state%th%active
435 u_qt_advection_enabled=is_component_field_available(
"u_advection") .and. &
436 is_component_field_available(
"th_advection") .and. us_qt_enabled
437 u_qt_viscosity_diffusion_enabled=is_component_field_available(
"u_viscosity") .and. &
438 is_component_field_available(
"th_diffusion") .and. us_qt_enabled
439 wu_qt_enabled=current_state%w%active .and. us_qt_enabled
443 if (us_qt_enabled)
then 445 allocate(
us_qt(column_size))
447 if (u_qt_advection_enabled)
then 451 if (u_qt_viscosity_diffusion_enabled)
then 455 if (wu_qt_enabled)
then 457 allocate(
wu_qt(column_size))
460 vs_qt_enabled=current_state%v%active .and. current_state%th%active
461 v_qt_advection_enabled=is_component_field_available(
"v_advection") .and. &
462 is_component_field_available(
"th_advection") .and. vs_qt_enabled
463 v_qt_viscosity_diffusion_enabled=is_component_field_available(
"v_viscosity") .and. &
464 is_component_field_available(
"th_diffusion") .and. vs_qt_enabled
465 wv_qt_enabled=current_state%w%active .and. vs_qt_enabled
469 if (vs_qt_enabled)
then 471 allocate(
vs_qt(column_size))
473 if (v_qt_advection_enabled)
then 477 if (v_qt_viscosity_diffusion_enabled)
then 481 if (wv_qt_enabled)
then 483 allocate(
wv_qt(column_size))
486 w_qt_enabled=current_state%w%active .and. current_state%th%active
487 ws_qt_enabled=w_qt_enabled
488 w_qt_advection_enabled=is_component_field_available(
"w_advection") .and. &
489 is_component_field_available(
"th_advection") .and. w_qt_enabled
490 w_qt_viscosity_diffusion_enabled=is_component_field_available(
"w_viscosity") .and. &
491 is_component_field_available(
"th_diffusion") .and. w_qt_enabled
492 w_qt_buoyancy_enabled=current_state%th%active .and. is_component_field_available(
"w_buoyancy")
493 ww_qt_enabled=w_qt_enabled
497 if (w_qt_enabled)
then 499 allocate(
w_qt(column_size))
501 if (ws_qt_enabled)
then 503 allocate(
ws_qt(column_size))
505 if (w_qt_advection_enabled)
then 509 if (w_qt_viscosity_diffusion_enabled)
then 513 if (w_qt_buoyancy_enabled)
then 517 if (ww_qt_enabled)
then 519 allocate(
ww_qt(column_size))
522 qt_qt_enabled=current_state%th%active
523 sqt_qt_enabled=qt_qt_enabled
524 qt_qt_advection_enabled=is_component_field_available(
"qt_advection") .and. qt_qt_enabled
525 qt_qt_diffusion_enabled=is_component_field_available(
"qt_diffusion") .and. qt_qt_enabled
526 wqt_qt_enabled=current_state%w%active .and. qt_qt_enabled
530 if (qt_qt_enabled)
then 532 allocate(
qt_qt(column_size))
534 if (sqt_qt_enabled)
then 536 allocate(
sqt_qt(column_size))
538 if (qt_qt_advection_enabled)
then 542 if (qt_qt_diffusion_enabled)
then 546 if (wqt_qt_enabled)
then 548 allocate(
wqt_qt(column_size))
556 type(model_state_type),
target,
intent(inout) :: current_state
558 integer :: column_size
559 logical :: u_mse_enabled, us_mse_enabled, u_mse_advection_enabled, u_mse_viscosity_diffusion_enabled, &
560 wu_mse_enabled, v_mse_enabled, vs_mse_enabled, v_mse_advection_enabled, &
561 v_mse_viscosity_diffusion_enabled, wv_mse_enabled, w_mse_enabled, ws_mse_enabled, &
562 w_mse_advection_enabled, w_mse_viscosity_diffusion_enabled, w_mse_buoyancy_enabled, ww_mse_enabled, &
563 mse_mse_enabled, smse_mse_enabled, mse_mse_advection_enabled, &
564 mse_mse_diffusion_enabled, wmse_mse_enabled
566 column_size=current_state%local_grid%size(z_index)
568 u_mse_enabled=current_state%u%active .and. current_state%th%active
569 us_mse_enabled=u_mse_enabled
570 u_mse_advection_enabled=is_component_field_available(
"u_advection") .and. &
571 is_component_field_available(
"th_advection") .and. u_mse_enabled
572 u_mse_viscosity_diffusion_enabled=is_component_field_available(
"u_viscosity") .and. &
573 is_component_field_available(
"th_diffusion") .and. u_mse_enabled
574 wu_mse_enabled=current_state%w%active .and. u_mse_enabled
578 if (u_mse_enabled)
then 580 allocate(
u_mse(column_size))
582 if (us_mse_enabled)
then 584 allocate(
us_mse(column_size))
586 if (u_mse_advection_enabled)
then 590 if (u_mse_viscosity_diffusion_enabled)
then 594 if (wu_mse_enabled)
then 596 allocate(
wu_mse(column_size))
599 v_mse_enabled=current_state%v%active .and. current_state%th%active
600 vs_mse_enabled=v_mse_enabled
601 v_mse_advection_enabled=is_component_field_available(
"v_advection") .and. &
602 is_component_field_available(
"th_advection") .and. v_mse_enabled
603 v_mse_viscosity_diffusion_enabled=is_component_field_available(
"v_viscosity") .and. &
604 is_component_field_available(
"th_diffusion") .and. v_mse_enabled
605 wv_mse_enabled=current_state%w%active .and. v_mse_enabled
609 if (v_mse_enabled)
then 611 allocate(
v_mse(column_size))
613 if (vs_mse_enabled)
then 615 allocate(
vs_mse(column_size))
617 if (v_mse_advection_enabled)
then 621 if (v_mse_viscosity_diffusion_enabled)
then 625 if (wv_mse_enabled)
then 627 allocate(
wv_mse(column_size))
630 w_mse_enabled=current_state%w%active .and. current_state%th%active
631 ws_mse_enabled=w_mse_enabled
632 w_mse_advection_enabled=is_component_field_available(
"w_advection") .and. &
633 is_component_field_available(
"th_advection") .and. w_mse_enabled
634 w_mse_viscosity_diffusion_enabled=is_component_field_available(
"w_viscosity") .and. &
635 is_component_field_available(
"th_diffusion") .and. w_mse_enabled
636 w_mse_buoyancy_enabled=current_state%th%active .and. is_component_field_available(
"w_buoyancy")
637 ww_mse_enabled=w_mse_enabled
641 if (w_mse_enabled)
then 643 allocate(
w_mse(column_size))
645 if (ws_mse_enabled)
then 647 allocate(
ws_mse(column_size))
649 if (w_mse_advection_enabled)
then 653 if (w_mse_viscosity_diffusion_enabled)
then 657 if (w_mse_buoyancy_enabled)
then 661 if (ww_mse_enabled)
then 663 allocate(
ww_mse(column_size))
666 mse_mse_enabled=current_state%th%active
667 smse_mse_enabled=mse_mse_enabled
668 mse_mse_advection_enabled=is_component_field_available(
"mse_advection") .and. mse_mse_enabled
669 mse_mse_diffusion_enabled=is_component_field_available(
"mse_diffusion") .and. mse_mse_enabled
670 wmse_mse_enabled=current_state%w%active .and. mse_mse_enabled
674 if (mse_mse_enabled)
then 678 if (smse_mse_enabled)
then 682 if (mse_mse_advection_enabled)
then 686 if (mse_mse_diffusion_enabled)
then 690 if (wmse_mse_enabled)
then 699 type(model_state_type),
target,
intent(inout) :: current_state
701 integer :: column_size
702 logical :: u_thetal_enabled, us_thetal_enabled, u_thetal_advection_enabled, u_thetal_viscosity_diffusion_enabled, &
703 wu_thetal_enabled, v_thetal_enabled, vs_thetal_enabled, v_thetal_advection_enabled, &
704 v_thetal_viscosity_diffusion_enabled, wv_thetal_enabled, w_thetal_enabled, ws_thetal_enabled, &
705 w_thetal_advection_enabled, w_thetal_viscosity_diffusion_enabled, w_thetal_buoyancy_enabled, ww_thetal_enabled, &
706 thetal_thetal_enabled, sthetal_thetal_enabled, thetal_thetal_advection_enabled, &
707 thetal_thetal_diffusion_enabled, wthetal_thetal_enabled
709 column_size=current_state%local_grid%size(z_index)
711 u_thetal_enabled=current_state%u%active .and. current_state%th%active
712 us_thetal_enabled=u_thetal_enabled
713 u_thetal_advection_enabled=is_component_field_available(
"u_advection") .and. &
714 is_component_field_available(
"th_advection") .and. u_thetal_enabled
715 u_thetal_viscosity_diffusion_enabled=is_component_field_available(
"u_viscosity") .and. &
716 is_component_field_available(
"th_diffusion") .and. u_thetal_enabled
717 wu_thetal_enabled=current_state%w%active .and. u_thetal_enabled
721 if (u_thetal_enabled)
then 725 if (us_thetal_enabled)
then 729 if (u_thetal_advection_enabled)
then 733 if (u_thetal_viscosity_diffusion_enabled)
then 737 if (wu_thetal_enabled)
then 742 v_thetal_enabled=current_state%v%active .and. current_state%th%active
743 vs_thetal_enabled=v_thetal_enabled
744 v_thetal_advection_enabled=is_component_field_available(
"v_advection") .and. &
745 is_component_field_available(
"th_advection") .and. v_thetal_enabled
746 v_thetal_viscosity_diffusion_enabled=is_component_field_available(
"v_viscosity") .and. &
747 is_component_field_available(
"th_diffusion") .and. v_thetal_enabled
748 wv_thetal_enabled=current_state%w%active .and. v_thetal_enabled
752 if (v_thetal_enabled)
then 756 if (vs_thetal_enabled)
then 760 if (v_thetal_advection_enabled)
then 764 if (v_thetal_viscosity_diffusion_enabled)
then 768 if (wv_thetal_enabled)
then 773 w_thetal_enabled=current_state%w%active .and. current_state%th%active
774 ws_thetal_enabled=w_thetal_enabled
775 w_thetal_advection_enabled=is_component_field_available(
"w_advection") .and. &
776 is_component_field_available(
"th_advection") .and. w_thetal_enabled
777 w_thetal_viscosity_diffusion_enabled=is_component_field_available(
"w_viscosity") .and. &
778 is_component_field_available(
"th_diffusion") .and. w_thetal_enabled
779 w_thetal_buoyancy_enabled=current_state%th%active .and. is_component_field_available(
"w_buoyancy")
780 ww_thetal_enabled=w_thetal_enabled
784 if (w_thetal_enabled)
then 788 if (ws_thetal_enabled)
then 792 if (w_thetal_advection_enabled)
then 796 if (w_thetal_viscosity_diffusion_enabled)
then 800 if (w_thetal_buoyancy_enabled)
then 804 if (ww_thetal_enabled)
then 809 thetal_thetal_enabled=current_state%th%active
810 sthetal_thetal_enabled=thetal_thetal_enabled
811 thetal_thetal_advection_enabled=is_component_field_available(
"th_advection") .and. thetal_thetal_enabled
812 thetal_thetal_diffusion_enabled=is_component_field_available(
"th_diffusion") .and. thetal_thetal_enabled
813 wthetal_thetal_enabled=current_state%w%active .and. thetal_thetal_enabled
817 if (thetal_thetal_enabled)
then 821 if (sthetal_thetal_enabled)
then 825 if (thetal_thetal_advection_enabled)
then 829 if (thetal_thetal_diffusion_enabled)
then 833 if (wthetal_thetal_enabled)
then 842 type(model_state_type),
target,
intent(inout) :: current_state
844 integer :: column_size
845 logical :: tu_su_enabled, uu_advection_enabled, uu_viscosity_enabled, wu_u_enabled, tv_sv_enabled, vv_advection_enabled, &
846 vv_viscosity_enabled, wv_v_enabled, tw_sw_enabled, ww_advection_enabled, ww_viscosity_enabled, ww_buoyancy_enabled
848 tu_su_enabled=current_state%u%active
849 uu_advection_enabled=is_component_field_available(
"u_advection") .and. current_state%u%active
850 uu_viscosity_enabled=is_component_field_available(
"u_viscosity") .and. current_state%u%active
851 wu_u_enabled=current_state%u%active .and. current_state%w%active
852 tv_sv_enabled=current_state%v%active
853 vv_advection_enabled=is_component_field_available(
"v_advection") .and. current_state%v%active
854 vv_viscosity_enabled=is_component_field_available(
"v_viscosity") .and. current_state%v%active
855 wv_v_enabled=current_state%v%active .and. current_state%w%active
856 tw_sw_enabled=current_state%w%active
857 ww_advection_enabled=is_component_field_available(
"w_advection") .and. current_state%w%active
858 ww_viscosity_enabled=is_component_field_available(
"w_viscosity") .and. current_state%w%active
859 ww_buoyancy_enabled=is_component_field_available(
"w_buoyancy") .and. current_state%w%active
863 column_size=current_state%local_grid%size(z_index)
865 if (tu_su_enabled)
then 867 allocate(
tu_su(column_size))
869 if (uu_advection_enabled)
then 873 if (uu_viscosity_enabled)
then 877 if (wu_u_enabled)
then 879 allocate(
wu_u(column_size))
881 if (tv_sv_enabled)
then 883 allocate(
tv_sv(column_size))
885 if (vv_advection_enabled)
then 889 if (vv_viscosity_enabled)
then 893 if (wv_v_enabled)
then 895 allocate(
wv_v(column_size))
897 if (tw_sw_enabled)
then 899 allocate(
tw_sw(column_size))
901 if (ww_advection_enabled)
then 905 if (ww_viscosity_enabled)
then 909 if (ww_buoyancy_enabled)
then 918 type(model_state_type),
target,
intent(inout) :: current_state
920 integer :: column_size
921 logical :: uw_advection_term_enabled, vw_advection_term_enabled, uw_viscosity_term_enabled, &
922 vw_viscosity_term_enabled, uw_buoyancy_term_enabled, vw_buoyancy_term_enabled, uw_tendency_term_enabled, &
923 vw_tendency_term_enabled, uw_w_term_enabled, vw_w_term_enabled
925 uw_advection_term_enabled=is_component_field_available(
"w_advection") .and. is_component_field_available(
"u_advection") &
926 .and. current_state%w%active
927 vw_advection_term_enabled=is_component_field_available(
"w_advection") .and. is_component_field_available(
"v_advection") &
928 .and. current_state%w%active
929 uw_viscosity_term_enabled=is_component_field_available(
"w_viscosity") .and. is_component_field_available(
"u_viscosity") &
930 .and. current_state%w%active
931 vw_viscosity_term_enabled=is_component_field_available(
"w_viscosity") .and. is_component_field_available(
"v_viscosity") &
932 .and. current_state%w%active
933 uw_buoyancy_term_enabled=is_component_field_available(
"w_buoyancy")
934 vw_buoyancy_term_enabled=is_component_field_available(
"w_buoyancy")
935 uw_tendency_term_enabled=current_state%w%active .and. current_state%u%active
936 vw_tendency_term_enabled=current_state%w%active .and. current_state%v%active
937 uw_w_term_enabled=current_state%w%active .and. current_state%u%active
938 vw_w_term_enabled=current_state%w%active .and. current_state%v%active
941 vw_w_term_enabled .or. uw_advection_term_enabled .or. vw_advection_term_enabled .or. uw_viscosity_term_enabled .or. &
942 vw_viscosity_term_enabled
944 column_size=current_state%local_grid%size(z_index)
946 if (uw_advection_term_enabled)
then 950 if (vw_advection_term_enabled)
then 954 if (uw_viscosity_term_enabled)
then 958 if (vw_viscosity_term_enabled)
then 962 if (uw_buoyancy_term_enabled)
then 966 if (vw_buoyancy_term_enabled)
then 970 if (uw_tendency_term_enabled)
then 974 if (vw_tendency_term_enabled)
then 978 if (uw_w_term_enabled)
then 980 allocate(
uw_w(column_size))
982 if (vw_w_term_enabled)
then 984 allocate(
vw_w(column_size))
991 type(model_state_type),
target,
intent(inout) :: current_state
993 logical :: q_flux_term_enabled, q_tendency_term_enabled, q_gradient_term_enabled, q_diff_enabled, &
996 q_flux_term_enabled=current_state%number_q_fields .gt. 0 .and. current_state%w%active
997 q_tendency_term_enabled=current_state%number_q_fields .gt. 0 .and. current_state%w%active
998 q_gradient_term_enabled=is_component_field_available(
"w_advection") .and. is_component_field_available(
"q_advection") &
999 .and. current_state%w%active .and. current_state%number_q_fields .gt. 0
1000 q_diff_enabled=is_component_field_available(
"q_diffusion") .and. is_component_field_available(
"w_viscosity") &
1001 .and. current_state%w%active .and. current_state%number_q_fields .gt. 0
1002 q_buoyancy_enabled=is_component_field_available(
"w_buoyancy") .and. current_state%number_q_fields .gt. 0
1006 if (q_flux_term_enabled)
then 1008 allocate(
q_flux_values(current_state%local_grid%size(z_index), current_state%number_q_fields))
1010 if (q_tendency_term_enabled)
then 1012 allocate(
q_tendency(current_state%local_grid%size(z_index), current_state%number_q_fields))
1014 if (q_gradient_term_enabled)
then 1016 allocate(
q_gradient(current_state%local_grid%size(z_index), current_state%number_q_fields))
1018 if (q_diff_enabled)
then 1020 allocate(
q_diff(current_state%local_grid%size(z_index), current_state%number_q_fields))
1022 if (q_buoyancy_enabled)
then 1024 allocate(
q_buoyancy(current_state%local_grid%size(z_index), current_state%number_q_fields))
1031 type(model_state_type),
target,
intent(inout) :: current_state
1033 logical :: th_flux_term_enabled, th_tendency_term_enabled, th_diff_enabled, th_gradient_term_enabled, th_buoyancy_enabled
1035 th_flux_term_enabled=current_state%th%active .and. current_state%w%active
1036 th_tendency_term_enabled=current_state%th%active .and. current_state%w%active
1037 th_gradient_term_enabled=is_component_field_available(
"w_advection") .and. is_component_field_available(
"th_advection") &
1038 .and. current_state%w%active .and. current_state%th%active
1039 th_diff_enabled=is_component_field_available(
"th_diffusion") .and. is_component_field_available(
"w_viscosity") &
1040 .and. current_state%w%active .and. current_state%th%active
1041 th_buoyancy_enabled=is_component_field_available(
"w_buoyancy") .and. current_state%th%active
1045 if (th_flux_term_enabled)
then 1049 if (th_tendency_term_enabled)
then 1051 allocate(
th_tendency(current_state%local_grid%size(z_index)))
1053 if (th_diff_enabled)
then 1055 allocate(
th_diff(current_state%local_grid%size(z_index)))
1057 if (th_gradient_term_enabled)
then 1059 allocate(
th_gradient(current_state%local_grid%size(z_index)))
1061 if (th_buoyancy_enabled)
then 1063 allocate(
th_buoyancy(current_state%local_grid%size(z_index)))
1072 if (
allocated(
wu_qt))
wu_qt=0.0_default_precision
1073 if (
allocated(
vs_qt))
vs_qt=0.0_default_precision
1076 if (
allocated(
wv_qt))
wv_qt=0.0_default_precision
1077 if (
allocated(
w_qt))
w_qt=0.0_default_precision
1078 if (
allocated(
ws_qt))
ws_qt=0.0_default_precision
1082 if (
allocated(
ww_qt))
ww_qt=0.0_default_precision
1083 if (
allocated(
qt_qt))
qt_qt=0.0_default_precision
1094 type(model_state_type),
target,
intent(inout) :: current_state
1096 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: upr, vpr, uprm1, vprm1, qtpr, qtprp1
1097 type(component_field_value_type) :: u_advection, u_viscosity, th_advection, th_diffusion, v_advection, v_viscosity, &
1098 w_advection, w_viscosity, w_buoyancy
1101 if (is_component_field_available(
"u_advection")) u_advection=get_component_field_value(current_state,
"u_advection")
1102 if (is_component_field_available(
"u_viscosity")) u_viscosity=get_component_field_value(current_state,
"u_viscosity")
1103 if (is_component_field_available(
"th_advection")) th_advection=get_component_field_value(current_state,
"th_advection")
1104 if (is_component_field_available(
"th_diffusion")) th_diffusion=get_component_field_value(current_state,
"th_diffusion")
1105 if (is_component_field_available(
"v_advection")) v_advection=get_component_field_value(current_state,
"v_advection")
1106 if (is_component_field_available(
"v_viscosity")) v_viscosity=get_component_field_value(current_state,
"v_viscosity")
1107 if (is_component_field_available(
"w_advection")) w_advection=get_component_field_value(current_state,
"w_advection")
1108 if (is_component_field_available(
"w_viscosity")) w_viscosity=get_component_field_value(current_state,
"w_viscosity")
1109 if (is_component_field_available(
"w_buoyancy")) w_buoyancy=get_component_field_value(current_state,
"w_buoyancy")
1111 do k=1, current_state%local_grid%size(z_index)
1112 upr(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x)
1113 uprm1(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x-1)
1114 if (
allocated(current_state%global_grid%configuration%vertical%olubar))
then 1115 upr(k)=upr(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1116 uprm1(k)=uprm1(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1118 vpr(k)=current_state%v%data(k,current_state%column_local_y,current_state%column_local_x)
1119 vprm1(k)=current_state%v%data(k,current_state%column_local_y-1,current_state%column_local_x)
1120 if (
allocated(current_state%global_grid%configuration%vertical%olvbar))
then 1121 vpr(k)=vpr(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1122 vprm1(k)=vprm1(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1125 qtpr(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x)
1126 qtprp1(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x+1)
1127 if (
allocated(current_state%global_grid%configuration%vertical%olthbar))
then 1128 qtpr(k)=qtpr(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1129 qtprp1(k)=qtprp1(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1132 do k=2, current_state%local_grid%size(z_index)-1
1134 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(qtpr(k)+qtprp1(k))*&
1135 current_state%su%data(k,current_state%column_local_y,current_state%column_local_x)
1137 th_advection%real_1d_array(k)+0.5*(qtpr(k)+qtprp1(k))*u_advection%real_1d_array(k)
1139 (qtpr(k)+qtprp1(k))*u_viscosity%real_1d_array(k)+0.5*(upr(k)+uprm1(k))*th_diffusion%real_1d_array(k)
1141 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.25*(upr(k+1)+upr(k)+&
1142 uprm1(k+1)+uprm1(k))*0.5*(qtpr(k+1)+qtpr(k))
1145 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(qtpr(k)+qtprp1(k))*&
1146 current_state%sv%data(k,current_state%column_local_y,current_state%column_local_x)
1148 th_advection%real_1d_array(k)+0.5*(qtpr(k)+qtprp1(k))*v_advection%real_1d_array(k)
1150 (qtpr(k)+qtprp1(k))*v_viscosity%real_1d_array(k)+0.5*(vpr(k)+vprm1(k))*th_diffusion%real_1d_array(k)
1152 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.25*(vpr(k+1)+vpr(k)+&
1153 vprm1(k+1)+vprm1(k))*0.5*(qtpr(k+1)+qtpr(k))
1156 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(qtpr(k)+qtpr(k+1))
1158 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*&
1159 (current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+&
1160 current_state%sth%data(k+1,current_state%column_local_y,current_state%column_local_x))+0.5*(qtpr(k)+qtpr(k+1))*&
1161 current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)
1163 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(th_advection%real_1d_array(k)+&
1164 th_advection%real_1d_array(k+1))+0.5*(qtpr(k)+qtpr(k+1))*w_advection%real_1d_array(k)
1166 (qtpr(k)+qtpr(k+1))*w_viscosity%real_1d_array(k)+&
1167 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*&
1168 (th_diffusion%real_1d_array(k)+th_diffusion%real_1d_array(k+1))
1170 w_buoyancy%real_1d_array(k)
1172 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1173 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*0.5*(&
1174 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1175 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*qtpr(k)
1179 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)
1181 th_advection%real_1d_array(k)
1183 th_diffusion%real_1d_array(k)
1185 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(qtpr(k+1)+qtpr(k))*0.5*(&
1189 if (
allocated(u_advection%real_1d_array))
deallocate(u_advection%real_1d_array)
1190 if (
allocated(u_viscosity%real_1d_array))
deallocate(u_viscosity%real_1d_array)
1191 if (
allocated(th_advection%real_1d_array))
deallocate(th_advection%real_1d_array)
1192 if (
allocated(th_diffusion%real_1d_array))
deallocate(th_diffusion%real_1d_array)
1193 if (
allocated(v_advection%real_1d_array))
deallocate(v_advection%real_1d_array)
1194 if (
allocated(v_viscosity%real_1d_array))
deallocate(v_viscosity%real_1d_array)
1195 if (
allocated(w_advection%real_1d_array))
deallocate(w_advection%real_1d_array)
1196 if (
allocated(w_viscosity%real_1d_array))
deallocate(w_viscosity%real_1d_array)
1197 if (
allocated(w_buoyancy%real_1d_array))
deallocate(w_buoyancy%real_1d_array)
1208 type(model_state_type),
target,
intent(inout) :: current_state
1212 do k=2, current_state%local_grid%size(z_index)-1
1213 if (current_state%w%data(k, current_state%column_local_y,current_state%column_local_x) .gt.
wmfcrit)
then 1214 mflux=
mflux+current_state%global_grid%configuration%vertical%rho(k)*&
1215 current_state%global_grid%configuration%vertical%dzn(k)*&
1216 current_state%w%data(k, current_state%column_local_y,current_state%column_local_x)
1228 if (
allocated(
v_mse))
v_mse=0.0_default_precision
1233 if (
allocated(
w_mse))
w_mse=0.0_default_precision
1250 type(model_state_type),
target,
intent(inout) :: current_state
1252 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: upr, vpr, uprm1, vprm1, msepr, mseprp1
1253 type(component_field_value_type) :: u_advection, u_viscosity, th_advection, th_diffusion, v_advection, v_viscosity, &
1254 w_advection, w_viscosity, w_buoyancy
1257 if (is_component_field_available(
"u_advection")) u_advection=get_component_field_value(current_state,
"u_advection")
1258 if (is_component_field_available(
"u_viscosity")) u_viscosity=get_component_field_value(current_state,
"u_viscosity")
1259 if (is_component_field_available(
"th_advection")) th_advection=get_component_field_value(current_state,
"th_advection")
1260 if (is_component_field_available(
"th_diffusion")) th_diffusion=get_component_field_value(current_state,
"th_diffusion")
1261 if (is_component_field_available(
"v_advection")) v_advection=get_component_field_value(current_state,
"v_advection")
1262 if (is_component_field_available(
"v_viscosity")) v_viscosity=get_component_field_value(current_state,
"v_viscosity")
1263 if (is_component_field_available(
"w_advection")) w_advection=get_component_field_value(current_state,
"w_advection")
1264 if (is_component_field_available(
"w_viscosity")) w_viscosity=get_component_field_value(current_state,
"w_viscosity")
1265 if (is_component_field_available(
"w_buoyancy")) w_buoyancy=get_component_field_value(current_state,
"w_buoyancy")
1267 do k=1, current_state%local_grid%size(z_index)
1268 upr(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x)
1269 uprm1(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x-1)
1270 if (
allocated(current_state%global_grid%configuration%vertical%olubar))
then 1271 upr(k)=upr(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1272 uprm1(k)=uprm1(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1274 vpr(k)=current_state%v%data(k,current_state%column_local_y,current_state%column_local_x)
1275 vprm1(k)=current_state%v%data(k,current_state%column_local_y-1,current_state%column_local_x)
1276 if (
allocated(current_state%global_grid%configuration%vertical%olvbar))
then 1277 vpr(k)=vpr(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1278 vprm1(k)=vprm1(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1280 msepr(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x)
1281 mseprp1(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x+1)
1282 if (
allocated(current_state%global_grid%configuration%vertical%olthbar))
then 1283 msepr(k)=msepr(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1284 mseprp1(k)=mseprp1(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1287 do k=2, current_state%local_grid%size(z_index)-1
1288 if (
allocated(
u_mse))
u_mse(k)=
u_mse(k)+0.5*(upr(k)+uprm1(k))*msepr(k)
1290 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(msepr(k)+mseprp1(k))*&
1291 current_state%su%data(k,current_state%column_local_y,current_state%column_local_x)
1293 th_advection%real_1d_array(k)+0.5*(msepr(k)+mseprp1(k))*u_advection%real_1d_array(k)
1295 (msepr(k)+mseprp1(k))*u_viscosity%real_1d_array(k)+0.5*(upr(k)+uprm1(k))*th_diffusion%real_1d_array(k)
1297 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.25*(upr(k+1)+upr(k)+&
1298 uprm1(k+1)+uprm1(k))*0.5*(msepr(k+1)+msepr(k))
1300 if (
allocated(
v_mse))
v_mse(k)=
v_mse(k)+0.5*(vpr(k)+vprm1(k))*msepr(k)
1302 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(msepr(k)+mseprp1(k))*&
1303 current_state%sv%data(k,current_state%column_local_y,current_state%column_local_x)
1305 th_advection%real_1d_array(k)+0.5*(msepr(k)+mseprp1(k))*v_advection%real_1d_array(k)
1307 (msepr(k)+mseprp1(k))*v_viscosity%real_1d_array(k)+0.5*(vpr(k)+vprm1(k))*th_diffusion%real_1d_array(k)
1309 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.25*(vpr(k+1)+vpr(k)+&
1310 vprm1(k+1)+vprm1(k))*0.5*(msepr(k+1)+msepr(k))
1313 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(msepr(k)+msepr(k+1))
1315 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*&
1316 (current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+&
1317 current_state%sth%data(k+1,current_state%column_local_y,current_state%column_local_x))+0.5*(msepr(k)+msepr(k+1))*&
1318 current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)
1320 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(th_advection%real_1d_array(k)+&
1321 th_advection%real_1d_array(k+1))+0.5*(msepr(k)+msepr(k+1))*w_advection%real_1d_array(k)
1323 (msepr(k)+msepr(k+1))*w_viscosity%real_1d_array(k)+&
1324 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*&
1325 (th_diffusion%real_1d_array(k)+th_diffusion%real_1d_array(k+1))
1327 w_buoyancy%real_1d_array(k)
1329 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1330 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*0.5*(&
1331 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1332 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*msepr(k)
1336 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)
1338 th_advection%real_1d_array(k)
1340 th_diffusion%real_1d_array(k)
1342 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(msepr(k+1)+msepr(k))*0.5*(&
1343 msepr(k+1)+msepr(k))
1346 if (
allocated(u_advection%real_1d_array))
deallocate(u_advection%real_1d_array)
1347 if (
allocated(u_viscosity%real_1d_array))
deallocate(u_viscosity%real_1d_array)
1348 if (
allocated(th_advection%real_1d_array))
deallocate(th_advection%real_1d_array)
1349 if (
allocated(th_diffusion%real_1d_array))
deallocate(th_diffusion%real_1d_array)
1350 if (
allocated(v_advection%real_1d_array))
deallocate(v_advection%real_1d_array)
1351 if (
allocated(v_viscosity%real_1d_array))
deallocate(v_viscosity%real_1d_array)
1352 if (
allocated(w_advection%real_1d_array))
deallocate(w_advection%real_1d_array)
1353 if (
allocated(w_viscosity%real_1d_array))
deallocate(w_viscosity%real_1d_array)
1354 if (
allocated(w_buoyancy%real_1d_array))
deallocate(w_buoyancy%real_1d_array)
1385 type(model_state_type),
target,
intent(inout) :: current_state
1387 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: upr, vpr, uprm1, vprm1, thlpr, thlprp1
1388 type(component_field_value_type) :: u_advection, u_viscosity, th_advection, th_diffusion, v_advection, v_viscosity, &
1389 w_advection, w_viscosity, w_buoyancy
1392 if (is_component_field_available(
"u_advection")) u_advection=get_component_field_value(current_state,
"u_advection")
1393 if (is_component_field_available(
"u_viscosity")) u_viscosity=get_component_field_value(current_state,
"u_viscosity")
1394 if (is_component_field_available(
"th_advection")) th_advection=get_component_field_value(current_state,
"th_advection")
1395 if (is_component_field_available(
"th_diffusion")) th_diffusion=get_component_field_value(current_state,
"th_diffusion")
1396 if (is_component_field_available(
"v_advection")) v_advection=get_component_field_value(current_state,
"v_advection")
1397 if (is_component_field_available(
"v_viscosity")) v_viscosity=get_component_field_value(current_state,
"v_viscosity")
1398 if (is_component_field_available(
"w_advection")) w_advection=get_component_field_value(current_state,
"w_advection")
1399 if (is_component_field_available(
"w_viscosity")) w_viscosity=get_component_field_value(current_state,
"w_viscosity")
1400 if (is_component_field_available(
"w_buoyancy")) w_buoyancy=get_component_field_value(current_state,
"w_buoyancy")
1402 do k=1, current_state%local_grid%size(z_index)
1403 upr(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x)
1404 uprm1(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x-1)
1405 if (
allocated(current_state%global_grid%configuration%vertical%olubar))
then 1406 upr(k)=upr(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1407 uprm1(k)=uprm1(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1409 vpr(k)=current_state%v%data(k,current_state%column_local_y,current_state%column_local_x)
1410 vprm1(k)=current_state%v%data(k,current_state%column_local_y-1,current_state%column_local_x)
1411 if (
allocated(current_state%global_grid%configuration%vertical%olvbar))
then 1412 vpr(k)=vpr(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1413 vprm1(k)=vprm1(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1415 thlpr(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x)
1416 thlprp1(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x+1)
1417 if (
allocated(current_state%global_grid%configuration%vertical%olthbar))
then 1418 thlpr(k)=thlpr(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1419 thlprp1(k)=thlprp1(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1422 do k=2, current_state%local_grid%size(z_index)-1
1425 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(thlpr(k)+thlprp1(k))*&
1426 current_state%su%data(k,current_state%column_local_y,current_state%column_local_x)
1428 th_advection%real_1d_array(k)+0.5*(thlpr(k)+thlprp1(k))*u_advection%real_1d_array(k)
1430 (thlpr(k)+thlprp1(k))*u_viscosity%real_1d_array(k)+0.5*(upr(k)+uprm1(k))*th_diffusion%real_1d_array(k)
1432 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.25*(upr(k+1)+upr(k)+&
1433 uprm1(k+1)+uprm1(k))*0.5*(thlpr(k+1)+thlpr(k))
1437 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(thlpr(k)+thlprp1(k))*&
1438 current_state%sv%data(k,current_state%column_local_y,current_state%column_local_x)
1440 th_advection%real_1d_array(k)+0.5*(thlpr(k)+thlprp1(k))*v_advection%real_1d_array(k)
1442 (thlpr(k)+thlprp1(k))*v_viscosity%real_1d_array(k)+0.5*(vpr(k)+vprm1(k))*th_diffusion%real_1d_array(k)
1444 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.25*(vpr(k+1)+vpr(k)+&
1445 vprm1(k+1)+vprm1(k))*0.5*(thlpr(k+1)+thlpr(k))
1448 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(thlpr(k)+thlpr(k+1))
1450 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*&
1451 (current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+&
1452 current_state%sth%data(k+1,current_state%column_local_y,current_state%column_local_x))+0.5*(thlpr(k)+thlpr(k+1))*&
1453 current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)
1455 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(th_advection%real_1d_array(k)+&
1456 th_advection%real_1d_array(k+1))+0.5*(thlpr(k)+thlpr(k+1))*w_advection%real_1d_array(k)
1458 (thlpr(k)+thlpr(k+1))*w_viscosity%real_1d_array(k)+&
1459 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*&
1460 (th_diffusion%real_1d_array(k)+th_diffusion%real_1d_array(k+1))
1462 w_buoyancy%real_1d_array(k)
1464 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1465 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*0.5*(&
1466 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1467 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*thlpr(k)
1471 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)
1473 th_advection%real_1d_array(k)
1475 th_diffusion%real_1d_array(k)
1477 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(thlpr(k+1)+thlpr(k))*0.5*(&
1478 thlpr(k+1)+thlpr(k))
1481 if (
allocated(u_advection%real_1d_array))
deallocate(u_advection%real_1d_array)
1482 if (
allocated(u_viscosity%real_1d_array))
deallocate(u_viscosity%real_1d_array)
1483 if (
allocated(th_advection%real_1d_array))
deallocate(th_advection%real_1d_array)
1484 if (
allocated(th_diffusion%real_1d_array))
deallocate(th_diffusion%real_1d_array)
1485 if (
allocated(v_advection%real_1d_array))
deallocate(v_advection%real_1d_array)
1486 if (
allocated(v_viscosity%real_1d_array))
deallocate(v_viscosity%real_1d_array)
1487 if (
allocated(w_advection%real_1d_array))
deallocate(w_advection%real_1d_array)
1488 if (
allocated(w_viscosity%real_1d_array))
deallocate(w_viscosity%real_1d_array)
1489 if (
allocated(w_buoyancy%real_1d_array))
deallocate(w_buoyancy%real_1d_array)
1497 if (
allocated(
wu_u))
wu_u=0.0_default_precision
1498 if (
allocated(
tv_sv))
tv_sv=0.0_default_precision
1501 if (
allocated(
wv_v))
wv_v=0.0_default_precision
1502 if (
allocated(
tw_sw))
tw_sw=0.0_default_precision
1511 type(model_state_type),
target,
intent(inout) :: current_state
1513 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: upr, vpr, uprm1, vprm1
1514 type(component_field_value_type) :: u_advection, u_viscosity, v_advection, v_viscosity, w_advection, w_viscosity, w_buoyancy
1517 if (is_component_field_available(
"u_advection")) u_advection=get_component_field_value(current_state,
"u_advection")
1518 if (is_component_field_available(
"u_viscosity")) u_viscosity=get_component_field_value(current_state,
"u_viscosity")
1519 if (is_component_field_available(
"v_advection")) v_advection=get_component_field_value(current_state,
"v_advection")
1520 if (is_component_field_available(
"v_viscosity")) v_viscosity=get_component_field_value(current_state,
"v_viscosity")
1521 if (is_component_field_available(
"w_advection")) w_advection=get_component_field_value(current_state,
"w_advection")
1522 if (is_component_field_available(
"w_viscosity")) w_viscosity=get_component_field_value(current_state,
"w_viscosity")
1523 if (is_component_field_available(
"w_buoyancy")) w_buoyancy=get_component_field_value(current_state,
"w_buoyancy")
1525 do k=1, current_state%local_grid%size(z_index)
1526 upr(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x)
1527 uprm1(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x-1)
1528 if (
allocated(current_state%global_grid%configuration%vertical%olubar))
then 1529 upr(k)=upr(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1530 uprm1(k)=uprm1(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1532 vpr(k)=current_state%v%data(k,current_state%column_local_y,current_state%column_local_x)
1533 vprm1(k)=current_state%v%data(k,current_state%column_local_y-1,current_state%column_local_x)
1534 if (
allocated(current_state%global_grid%configuration%vertical%olvbar))
then 1535 vpr(k)=vpr(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1536 vprm1(k)=vprm1(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1539 do k=2, current_state%local_grid%size(z_index)-1
1540 if (
allocated(
tu_su))
tu_su(k)=
tu_su(k)+2.0*upr(k)*current_state%su%data(k,current_state%column_local_y,&
1541 current_state%column_local_x)
1544 if (
allocated(
wu_u))
wu_u(k)=
wu_u(k)+0.25*(upr(k)+upr(k+1)+uprm1(k)+uprm1(k+1))*0.25*&
1545 (upr(k)+upr(k+1)+uprm1(k)+uprm1(k+1))*current_state%w%data(k,current_state%column_local_y,&
1546 current_state%column_local_x)
1547 if (
allocated(
tv_sv))
tv_sv(k)=
tv_sv(k)+2.0*vpr(k)*current_state%sv%data(k,current_state%column_local_y,&
1548 current_state%column_local_x)
1551 if (
allocated(
wv_v))
wv_v(k)=
wv_v(k)+0.25*(vpr(k)+vpr(k+1)+vprm1(k)+vprm1(k+1))*0.25*&
1552 (vpr(k)+vpr(k+1)+vprm1(k)+vprm1(k+1))*current_state%w%data(k,current_state%column_local_y,&
1553 current_state%column_local_x)
1554 if (
allocated(
tw_sw))
tw_sw(k)=
tw_sw(k)+2.0*current_state%w%data(k,current_state%column_local_y,&
1555 current_state%column_local_x)*current_state%sw%data(k,current_state%column_local_y,&
1556 current_state%column_local_x)
1558 current_state%column_local_x)*w_advection%real_1d_array(k)
1560 current_state%column_local_x)*w_viscosity%real_1d_array(k)
1562 current_state%column_local_x)*w_buoyancy%real_1d_array(k)
1565 if (
allocated(u_advection%real_1d_array))
deallocate(u_advection%real_1d_array)
1566 if (
allocated(u_viscosity%real_1d_array))
deallocate(u_viscosity%real_1d_array)
1567 if (
allocated(v_advection%real_1d_array))
deallocate(v_advection%real_1d_array)
1568 if (
allocated(v_viscosity%real_1d_array))
deallocate(v_viscosity%real_1d_array)
1569 if (
allocated(w_advection%real_1d_array))
deallocate(w_advection%real_1d_array)
1570 if (
allocated(w_viscosity%real_1d_array))
deallocate(w_viscosity%real_1d_array)
1571 if (
allocated(w_buoyancy%real_1d_array))
deallocate(w_buoyancy%real_1d_array)
1585 if (
allocated(
uw_w))
uw_w=0.0_default_precision
1586 if (
allocated(
vw_w))
vw_w=0.0_default_precision
1592 type(model_state_type),
target,
intent(inout) :: current_state
1594 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: upr, vpr, uprm1, vprm1
1595 type(component_field_value_type) :: w_advection, v_advection, u_advection, w_viscosity, v_viscosity, u_viscosity, &
1599 if (is_component_field_available(
"w_advection")) w_advection=get_component_field_value(current_state,
"w_advection")
1600 if (is_component_field_available(
"v_advection")) v_advection=get_component_field_value(current_state,
"v_advection")
1601 if (is_component_field_available(
"u_advection")) u_advection=get_component_field_value(current_state,
"u_advection")
1602 if (is_component_field_available(
"w_viscosity")) w_viscosity=get_component_field_value(current_state,
"w_viscosity")
1603 if (is_component_field_available(
"v_viscosity")) v_viscosity=get_component_field_value(current_state,
"v_viscosity")
1604 if (is_component_field_available(
"u_viscosity")) u_viscosity=get_component_field_value(current_state,
"u_viscosity")
1605 if (is_component_field_available(
"w_buoyancy")) w_buoyancy=get_component_field_value(current_state,
"w_buoyancy")
1607 do k=1, current_state%local_grid%size(z_index)
1608 upr(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x)
1609 uprm1(k)=current_state%u%data(k,current_state%column_local_y,current_state%column_local_x-1)
1610 if (
allocated(current_state%global_grid%configuration%vertical%olubar))
then 1611 upr(k)=upr(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1612 uprm1(k)=uprm1(k)-(current_state%global_grid%configuration%vertical%olubar(k)-current_state%ugal)
1614 vpr(k)=current_state%v%data(k,current_state%column_local_y,current_state%column_local_x)
1615 vprm1(k)=current_state%v%data(k,current_state%column_local_y-1,current_state%column_local_x)
1616 if (
allocated(current_state%global_grid%configuration%vertical%olvbar))
then 1617 vpr(k)=vpr(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1618 vprm1(k)=vprm1(k)-(current_state%global_grid%configuration%vertical%olvbar(k)-current_state%vgal)
1621 do k=2, current_state%local_grid%size(z_index)-1
1623 (w_advection%real_1d_array(k)+w_advection%real_1d_array(k-1))+0.25*(&
1624 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1625 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x)+&
1626 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x+1)+&
1627 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x+1))*u_advection%real_1d_array(k)
1629 (w_advection%real_1d_array(k)+w_advection%real_1d_array(k-1))+0.25*(&
1630 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1631 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x)+&
1632 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x+1)+&
1633 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x+1))*v_advection%real_1d_array(k)
1635 (w_viscosity%real_1d_array(k)+w_viscosity%real_1d_array(k-1))+0.25*(&
1636 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1637 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x)+&
1638 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x+1)+&
1639 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x+1))*u_viscosity%real_1d_array(k)
1641 (w_viscosity%real_1d_array(k)+w_viscosity%real_1d_array(k-1))+0.25*(&
1642 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1643 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x)+&
1644 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x+1)+&
1645 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x+1))*v_viscosity%real_1d_array(k)
1647 w_buoyancy%real_1d_array(k)+w_buoyancy%real_1d_array(k-1))
1649 w_buoyancy%real_1d_array(k)+w_buoyancy%real_1d_array(k-1))
1651 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1652 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x)+&
1653 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x+1)+&
1654 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x+1))*&
1655 current_state%su%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(&
1656 current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)+&
1657 current_state%sw%data(k-1,current_state%column_local_y,current_state%column_local_x))*0.5*(upr(k)+uprm1(k))
1659 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1660 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x)+&
1661 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x+1)+&
1662 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x+1))*&
1663 current_state%sv%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*(&
1664 current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)+&
1665 current_state%sw%data(k-1,current_state%column_local_y,current_state%column_local_x))*0.5*(vpr(k)+vprm1(k))
1666 if (
allocated(
uw_w))
uw_w(k)=
uw_w(k)+0.25*(upr(k)+upr(k+1)+uprm1(k)+uprm1(k+1))*&
1667 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*&
1668 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)
1669 if (
allocated(
vw_w))
vw_w(k)=
vw_w(k)+0.25*(vpr(k)+vpr(k+1)+vprm1(k)+vprm1(k+1))*&
1670 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*&
1671 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)
1674 if (
allocated(u_advection%real_1d_array))
deallocate(u_advection%real_1d_array)
1675 if (
allocated(u_viscosity%real_1d_array))
deallocate(u_viscosity%real_1d_array)
1676 if (
allocated(v_advection%real_1d_array))
deallocate(v_advection%real_1d_array)
1677 if (
allocated(v_viscosity%real_1d_array))
deallocate(v_viscosity%real_1d_array)
1678 if (
allocated(w_advection%real_1d_array))
deallocate(w_advection%real_1d_array)
1679 if (
allocated(w_viscosity%real_1d_array))
deallocate(w_viscosity%real_1d_array)
1680 if (
allocated(w_buoyancy%real_1d_array))
deallocate(w_buoyancy%real_1d_array)
1695 type(model_state_type),
target,
intent(inout) :: current_state
1698 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: qpr
1699 type(component_field_value_type) :: w_advection_published_value, q_advection_published_value, w_viscosity_published_value, &
1700 q_diffusion_published_value, w_buoyancy_published_value
1703 w_advection_published_value=get_component_field_value(current_state,
"w_advection")
1704 q_advection_published_value=get_component_field_value(current_state,
"q_advection")
1706 if (
allocated(
q_diff))
then 1707 w_viscosity_published_value=get_component_field_value(current_state,
"w_viscosity")
1708 q_diffusion_published_value=get_component_field_value(current_state,
"q_diffusion")
1711 w_buoyancy_published_value=get_component_field_value(current_state,
"w_buoyancy")
1714 do n=1, current_state%number_q_fields
1715 do k=1, current_state%local_grid%size(z_index)
1716 qpr(k)=current_state%q(n)%data(k,current_state%column_local_y,current_state%column_local_x)
1717 if (
allocated(current_state%global_grid%configuration%vertical%olqbar))
then 1718 qpr(k)=qpr(k)-current_state%global_grid%configuration%vertical%olqbar(k,n)
1721 do k=2, current_state%local_grid%size(z_index)-1
1723 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*&
1724 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(qpr(k)+qpr(k+1))
1726 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1727 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*&
1728 current_state%sq(n)%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*&
1729 qpr(k)*(current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)+&
1730 current_state%sw%data(k-1,current_state%column_local_y,current_state%column_local_x))
1733 w_advection_published_value%real_1d_array(k-1))+0.5*&
1734 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1735 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*&
1736 q_advection_published_value%real_2d_array(k, n)
1738 if (
allocated(
q_diff))
then 1739 q_diff(k,n)=
q_diff(k,n)+qpr(k)*0.5*(w_viscosity_published_value%real_1d_array(k)+&
1740 w_viscosity_published_value%real_1d_array(k-1))+&
1741 0.5*(current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1742 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*&
1743 q_diffusion_published_value%real_2d_array(k,n)
1747 w_buoyancy_published_value%real_1d_array(k-1))
1751 if (
allocated(w_advection_published_value%real_1d_array))
deallocate(w_advection_published_value%real_1d_array)
1752 if (
allocated(q_advection_published_value%real_2d_array))
deallocate(q_advection_published_value%real_2d_array)
1753 if (
allocated(w_viscosity_published_value%real_1d_array))
deallocate(w_viscosity_published_value%real_1d_array)
1754 if (
allocated(q_diffusion_published_value%real_2d_array))
deallocate(q_diffusion_published_value%real_2d_array)
1755 if (
allocated(w_buoyancy_published_value%real_1d_array))
deallocate(w_buoyancy_published_value%real_1d_array)
1770 type(model_state_type),
target,
intent(inout) :: current_state
1773 real(kind=DEFAULT_PRECISION),
dimension(current_state%local_grid%size(Z_INDEX)) :: thpr
1774 type(component_field_value_type) :: w_advection_published_value, th_advection_published_value, w_viscosity_published_value, &
1775 th_diffusion_published_value, w_buoyancy_published_value
1777 do k=1, current_state%local_grid%size(z_index)
1778 thpr(k)=current_state%th%data(k,current_state%column_local_y,current_state%column_local_x)
1779 if (
allocated(current_state%global_grid%configuration%vertical%olthbar))
then 1780 thpr(k)=thpr(k)-current_state%global_grid%configuration%vertical%olthbar(k)
1784 w_advection_published_value=get_component_field_value(current_state,
"w_advection")
1785 th_advection_published_value=get_component_field_value(current_state,
"th_advection")
1788 w_viscosity_published_value=get_component_field_value(current_state,
"w_viscosity")
1789 th_diffusion_published_value=get_component_field_value(current_state,
"th_diffusion")
1792 w_buoyancy_published_value=get_component_field_value(current_state,
"w_buoyancy")
1794 do k=2, current_state%local_grid%size(z_index)-1
1796 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*&
1797 current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)*0.5*(thpr(k)+thpr(k+1))
1799 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1800 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*&
1801 current_state%sth%data(k,current_state%column_local_y,current_state%column_local_x)+0.5*&
1802 thpr(k)*(current_state%sw%data(k,current_state%column_local_y,current_state%column_local_x)+&
1803 current_state%sw%data(k-1,current_state%column_local_y,current_state%column_local_x))
1806 w_advection_published_value%real_1d_array(k-1))+0.5*&
1807 (current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1808 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*&
1809 th_advection_published_value%real_1d_array(k)
1812 th_diff(k)=
th_diff(k)+thpr(k)*0.5*(w_viscosity_published_value%real_1d_array(k)+&
1813 w_viscosity_published_value%real_1d_array(k-1))+&
1814 0.5*(current_state%w%data(k,current_state%column_local_y,current_state%column_local_x)+&
1815 current_state%w%data(k-1,current_state%column_local_y,current_state%column_local_x))*&
1816 th_diffusion_published_value%real_1d_array(k)
1820 w_buoyancy_published_value%real_1d_array(k-1))
1823 if (
allocated(w_advection_published_value%real_1d_array))
deallocate(w_advection_published_value%real_1d_array)
1824 if (
allocated(th_advection_published_value%real_1d_array))
deallocate(th_advection_published_value%real_1d_array)
1825 if (
allocated(w_viscosity_published_value%real_1d_array))
deallocate(w_viscosity_published_value%real_1d_array)
1826 if (
allocated(th_diffusion_published_value%real_1d_array))
deallocate(th_diffusion_published_value%real_1d_array)
1827 if (
allocated(w_buoyancy_published_value%real_1d_array))
deallocate(w_buoyancy_published_value%real_1d_array)
1834 character(len=*),
intent(in) :: name
1843 character(len=*),
intent(in) :: name
1852 character(len=*),
intent(in) :: name
1861 character(len=*),
intent(in) :: name
1870 character(len=*),
intent(in) :: name
1879 character(len=*),
intent(in) :: name
1888 character(len=*),
intent(in) :: name
1897 character(len=*),
intent(in) :: name
1907 type(model_state_type),
target,
intent(inout) :: current_state
1908 character(len=*),
intent(in) :: name
1909 type(component_field_value_type),
intent(out) :: field_value
1911 if (name .eq.
"heat_flux_transport_local" .and.
allocated(
th_flux_values))
then 1913 else if (name .eq.
"heat_flux_gradient_local" .and.
allocated(
th_gradient))
then 1915 else if (name .eq.
"heat_flux_dissipation_local" .and.
allocated(
th_diff))
then 1917 else if (name .eq.
"heat_flux_buoyancy_local" .and.
allocated(
th_buoyancy))
then 1919 else if (name .eq.
"heat_flux_tendency_local" .and.
allocated(
th_tendency))
then 1921 else if (name .eq.
"q_flux_transport_local" .and.
allocated(
q_flux_values))
then 1923 else if (name .eq.
"q_flux_gradient_local" .and.
allocated(
q_gradient))
then 1925 else if (name .eq.
"q_flux_dissipation_local" .and.
allocated(
q_diff))
then 1927 else if (name .eq.
"q_flux_buoyancy_local" .and.
allocated(
q_buoyancy))
then 1929 else if (name .eq.
"q_flux_tendency_local" .and.
allocated(
q_tendency))
then 1931 else if (name .eq.
"uw_advection_local" .and.
allocated(
uw_advection))
then 1933 else if (name .eq.
"vw_advection_local" .and.
allocated(
vw_advection))
then 1935 else if (name .eq.
"uw_viscosity_local" .and.
allocated(
uw_viscosity))
then 1937 else if (name .eq.
"vw_viscosity_local" .and.
allocated(
vw_viscosity))
then 1939 else if (name .eq.
"uw_buoyancy_local" .and.
allocated(
uw_buoyancy))
then 1941 else if (name .eq.
"vw_buoyancy_local" .and.
allocated(
vw_buoyancy))
then 1943 else if (name .eq.
"uw_tendency_local" .and.
allocated(
uw_tendency))
then 1945 else if (name .eq.
"vw_tendency_local" .and.
allocated(
vw_tendency))
then 1947 else if (name .eq.
"uw_w_local" .and.
allocated(
uw_w))
then 1949 else if (name .eq.
"vw_w_local" .and.
allocated(
vw_w))
then 1951 else if (name .eq.
"tu_su_local" .and.
allocated(
tu_su))
then 1953 else if (name .eq.
"uu_advection_local" .and.
allocated(
uu_advection))
then 1955 else if (name .eq.
"uu_viscosity_local" .and.
allocated(
uu_viscosity))
then 1957 else if (name .eq.
"wu_u_local" .and.
allocated(
wu_u))
then 1959 else if (name .eq.
"tv_sv_local" .and.
allocated(
tv_sv))
then 1961 else if (name .eq.
"vv_advection_local" .and.
allocated(
vv_advection))
then 1963 else if (name .eq.
"vv_viscosity_local" .and.
allocated(
vv_viscosity))
then 1965 else if (name .eq.
"wv_v_local" .and.
allocated(
wv_v))
then 1967 else if (name .eq.
"tw_sw_local" .and.
allocated(
tw_sw))
then 1969 else if (name .eq.
"ww_advection_local" .and.
allocated(
ww_advection))
then 1971 else if (name .eq.
"ww_viscosity_local" .and.
allocated(
ww_viscosity))
then 1973 else if (name .eq.
"ww_buoyancy_local" .and.
allocated(
ww_buoyancy))
then 1975 else if (name .eq.
"u_thetal_local" .and.
allocated(
u_thetal))
then 1977 else if (name .eq.
"us_thetal_local" .and.
allocated(
us_thetal))
then 1979 else if (name .eq.
"u_thetal_advection_local" .and.
allocated(
u_thetal_advection))
then 1983 else if (name .eq.
"wu_thetal_local" .and.
allocated(
wu_thetal))
then 1985 else if (name .eq.
"v_thetal_local" .and.
allocated(
v_thetal))
then 1987 else if (name .eq.
"vs_thetal_local" .and.
allocated(
vs_thetal))
then 1989 else if (name .eq.
"v_thetal_advection_local" .and.
allocated(
v_thetal_advection))
then 1993 else if (name .eq.
"wv_thetal_local" .and.
allocated(
wv_thetal))
then 1995 else if (name .eq.
"w_thetal_local" .and.
allocated(
w_thetal))
then 1997 else if (name .eq.
"ws_thetal_local" .and.
allocated(
ws_thetal))
then 1999 else if (name .eq.
"w_thetal_advection_local" .and.
allocated(
w_thetal_advection))
then 2003 else if (name .eq.
"w_thetal_buoyancy_local" .and.
allocated(
w_thetal_buoyancy))
then 2005 else if (name .eq.
"ww_thetal_local" .and.
allocated(
ww_thetal))
then 2007 else if (name .eq.
"thetal_thetal_local" .and.
allocated(
thetal_thetal))
then 2009 else if (name .eq.
"sthetal_thetal_local" .and.
allocated(
sthetal_thetal))
then 2015 else if (name .eq.
"wthetal_thetal_local" .and.
allocated(
wthetal_thetal))
then 2017 else if (name .eq.
"u_mse_local" .and.
allocated(
u_mse))
then 2019 else if (name .eq.
"us_mse_local" .and.
allocated(
us_mse))
then 2021 else if (name .eq.
"u_mse_advection_local" .and.
allocated(
u_mse_advection))
then 2025 else if (name .eq.
"wu_mse_local" .and.
allocated(
wu_mse))
then 2027 else if (name .eq.
"v_mse_local" .and.
allocated(
v_mse))
then 2029 else if (name .eq.
"vs_mse_local" .and.
allocated(
vs_mse))
then 2031 else if (name .eq.
"v_mse_advection_local" .and.
allocated(
v_mse_advection))
then 2035 else if (name .eq.
"wv_mse_local" .and.
allocated(
wv_mse))
then 2037 else if (name .eq.
"w_mse_local" .and.
allocated(
w_mse))
then 2039 else if (name .eq.
"ws_mse_local" .and.
allocated(
ws_mse))
then 2041 else if (name .eq.
"w_mse_advection_local" .and.
allocated(
w_mse_advection))
then 2045 else if (name .eq.
"w_mse_buoyancy_local" .and.
allocated(
w_mse_buoyancy))
then 2047 else if (name .eq.
"ww_mse_local" .and.
allocated(
ww_mse))
then 2049 else if (name .eq.
"mse_mse_local" .and.
allocated(
mse_mse))
then 2051 else if (name .eq.
"smse_mse_local" .and.
allocated(
smse_mse))
then 2053 else if (name .eq.
"mse_mse_advection_local" .and.
allocated(
mse_mse_advection))
then 2055 else if (name .eq.
"mse_mse_diffusion_local" .and.
allocated(
mse_mse_diffusion))
then 2057 else if (name .eq.
"wmse_mse_local" .and.
allocated(
wmse_mse))
then 2059 else if (name .eq.
"us_qt_local" .and.
allocated(
us_qt))
then 2061 else if (name .eq.
"u_qt_advection_local" .and.
allocated(
u_qt_advection))
then 2065 else if (name .eq.
"wu_qt_local" .and.
allocated(
wu_qt))
then 2067 else if (name .eq.
"vs_qt_local" .and.
allocated(
vs_qt))
then 2069 else if (name .eq.
"v_qt_advection_local" .and.
allocated(
v_qt_advection))
then 2073 else if (name .eq.
"wv_qt_local" .and.
allocated(
wv_qt))
then 2075 else if (name .eq.
"w_qt_local" .and.
allocated(
w_qt))
then 2077 else if (name .eq.
"ws_qt_local" .and.
allocated(
ws_qt))
then 2079 else if (name .eq.
"w_qt_advection_local" .and.
allocated(
w_qt_advection))
then 2083 else if (name .eq.
"w_qt_buoyancy_local" .and.
allocated(
w_qt_buoyancy))
then 2085 else if (name .eq.
"ww_qt_local" .and.
allocated(
ww_qt))
then 2087 else if (name .eq.
"qt_qt_local" .and.
allocated(
qt_qt))
then 2089 else if (name .eq.
"sqt_qt_local" .and.
allocated(
sqt_qt))
then 2091 else if (name .eq.
"qt_qt_advection_local" .and.
allocated(
qt_qt_advection))
then 2093 else if (name .eq.
"qt_qt_diffusion_local" .and.
allocated(
qt_qt_diffusion))
then 2095 else if (name .eq.
"wqt_qt_local" .and.
allocated(
wqt_qt))
then 2097 else if (name .eq.
"mflux_local")
then 2098 field_value%scalar_real=
mflux 2107 type(component_field_value_type),
intent(inout) :: field_value
2108 real(kind=DEFAULT_PRECISION),
dimension(:),
optional :: real_1d_field
2109 real(kind=DEFAULT_PRECISION),
dimension(:,:),
optional :: real_2d_field
2111 if (
present(real_1d_field))
then 2112 allocate(field_value%real_1d_array(
size(real_1d_field)), source=real_1d_field)
2113 else if (
present(real_2d_field))
then 2114 allocate(field_value%real_2d_array(
size(real_2d_field, 1),
size(real_2d_field, 2)), source=real_2d_field)
2123 type(hashmap_type),
intent(inout) :: collection
2124 character(len=*),
intent(in) :: field_name
2125 logical,
intent(in) :: enabled_state
2127 call c_put_logical(collection, field_name, enabled_state)
2135 type(hashmap_type),
intent(inout) :: collection
2136 character(len=*),
intent(in) :: field_name
real(kind=default_precision), dimension(:), allocatable ww_qt
Gets a specific logical element out of the list, stack, queue or map with the corresponding key...
subroutine initialise_thetal_diagnostics(current_state)
Initialises the thetal diagnostics.
logical some_uw_vw_diagnostics_enabled
logical some_mse_diagnostics_enabled
type(hashmap_type) q_flux_fields
real(kind=default_precision), dimension(:), allocatable v_mse_advection
real(kind=default_precision), dimension(:), allocatable wv_qt
real(kind=default_precision), dimension(:), allocatable u_thetal
real(kind=default_precision), dimension(:), allocatable mse_mse_advection
real(kind=default_precision), dimension(:), allocatable tu_su
type(hashmap_type) prognostic_budget_fields
Wrapper type for the value returned for a published field from a component.
real(kind=default_precision), dimension(:), allocatable th_gradient
real(kind=default_precision), dimension(:), allocatable ww_buoyancy
real(kind=default_precision), dimension(:), allocatable u_mse_viscosity_diffusion
real(kind=default_precision), dimension(:), allocatable u_qt_advection
real(kind=default_precision), dimension(:), allocatable mse_mse
real(kind=default_precision), dimension(:), allocatable th_flux_values
real(kind=default_precision), dimension(:), allocatable u_thetal_viscosity_diffusion
real(kind=default_precision), dimension(:), allocatable uw_advection
subroutine clear_thetal()
Clears the thetal diagnostics.
real(kind=default_precision), dimension(:), allocatable u_qt_viscosity_diffusion
real(kind=default_precision) wmfcrit
real(kind=default_precision), dimension(:), allocatable vw_viscosity
real(kind=default_precision), dimension(:), allocatable th_diff
subroutine compute_q_flux_for_column(current_state)
Computes the Q flux diagnostics for a specific column.
logical function is_field_qt(name)
Determines whether a specific published field is a mse field.
real(kind=default_precision), dimension(:), allocatable u_mse
subroutine clear_uw_vw()
Clears the uw uv diagnostics.
real(kind=default_precision), dimension(:), allocatable w_qt_advection
subroutine set_published_field_value(field_value, real_1d_field, real_2d_field)
Sets the published field value from the temporary diagnostic values held by this component.
type(hashmap_type) mse_fields
real(kind=default_precision), dimension(:), allocatable ww_mse
real(kind=default_precision), dimension(:), allocatable us_mse
logical some_thetal_diagnostics_enabled
real(kind=default_precision), dimension(:), allocatable ws_thetal
subroutine field_information_retrieval_callback(current_state, name, field_information)
Field information retrieval callback, this returns information for a specific components published fi...
subroutine clear_qt()
Clears the qt diagnostics.
type(hashmap_type) qt_fields
real(kind=default_precision), dimension(:), allocatable uu_viscosity
logical function, public is_component_field_available(name)
Determines whether a specific published field is available or not.
subroutine initialise_mse_diagnostics(current_state)
Initialises the mse diagnostics. For now we are assuming mse is the same as theta, which needs updating with moisture information.
integer, parameter, public default_precision
MPI communication type which we use for the prognostic and calculation data.
integer, parameter, public z_index
Grid index parameters.
subroutine initialisation_callback(current_state)
Initialisation call back.
real(kind=default_precision), dimension(:), allocatable v_thetal_viscosity_diffusion
subroutine set_published_field_enabled_state(collection, field_name, enabled_state)
Sets the published value enabled state in the provided collection map.
real(kind=default_precision), dimension(:), allocatable w_thetal_viscosity_diffusion
real(kind=default_precision), dimension(:), allocatable vs_thetal
Contains common definitions for the data and datatypes used by MONC.
subroutine compute_scalars_for_column(current_state)
Computes the scalar diagnostics for a specific column.
real(kind=default_precision), dimension(:), allocatable v_mse_viscosity_diffusion
The ModelState which represents the current state of a run.
real(kind=default_precision), dimension(:), allocatable w_mse
real(kind=default_precision), dimension(:), allocatable vv_viscosity
real(kind=default_precision), dimension(:), allocatable wu_thetal
A hashmap structure, the same as a map but uses hashing for greatly improved performance when storing...
real(kind=default_precision), dimension(:), allocatable thetal_thetal_advection
real(kind=default_precision), dimension(:), allocatable w_qt
real(kind=default_precision), dimension(:), allocatable uw_w
subroutine initialise_scalar_diagnostics(current_state)
Initialises the scalar diagnostics.
real(kind=default_precision), dimension(:), allocatable wu_qt
subroutine initialise_qt_diagnostics(current_state)
Initialises the qt diagnostics. For now we are assuming qt is the same as theta, which needs updating...
Description of a component.
logical function is_field_scalar(name)
Determines whether a specific published field is a scalar field.
real(kind=default_precision), dimension(:), allocatable smse_mse
real(kind=default_precision), dimension(:), allocatable v_thetal_advection
logical some_theta_flux_diagnostics_enabled
subroutine clear_prognostic_budgets()
Clears the prognostic (uu, vv, ww) budgets.
real(kind=default_precision), dimension(:), allocatable vv_advection
subroutine initialise_uw_vw_diagnostics(current_state)
Initialises the UW and VW diagnostics.
real(kind=default_precision), dimension(:), allocatable ws_mse
real(kind=default_precision), dimension(:), allocatable us_thetal
real(kind=default_precision), dimension(:), allocatable w_thetal
integer, parameter, public component_array_field_type
real(kind=default_precision), dimension(:), allocatable vw_advection
real(kind=default_precision), dimension(:), allocatable qt_qt_advection
real(kind=default_precision), dimension(:), allocatable w_mse_buoyancy
Puts a logical key-value pair into the map.
subroutine clear_theta_fluxes()
Clears the heat flux diagnostics at the start of a timestep.
type(component_field_value_type) function, public get_component_field_value(current_state, name)
Retrieves the value wrapper of a components published field.
real(kind=default_precision), dimension(:), allocatable w_thetal_buoyancy
real(kind=default_precision), dimension(:), allocatable w_mse_viscosity_diffusion
subroutine initialise_q_flux_diagnostics(current_state)
Initialises the Q field flux diagnostic areas and enabled flags depending upon the configuration of t...
subroutine clear_scalars()
Clears the scalar diagnostics.
real(kind=default_precision), dimension(:,:), allocatable q_tendency
subroutine initialise_theta_flux_diagnostics(current_state)
Initialises the heat flux diagnostic areas and enabled flags depending upon the configuration of the ...
subroutine compute_thetal_for_column(current_state)
Computes the thetal diagnostics for a specific column.
real(kind=default_precision), dimension(:), allocatable us_qt
type(component_descriptor_type) function, public flux_budget_get_descriptor()
Provides the descriptor back to the caller and is used in component registration. ...
type(hashmap_type) thetal_fields
Returns the number of elements in the collection.
real(kind=default_precision), dimension(:), allocatable v_qt_viscosity_diffusion
real(kind=default_precision), dimension(:), allocatable w_qt_viscosity_diffusion
type(component_field_information_type) function, public get_component_field_information(current_state, name)
Retrieves information about a components published field which includes its type and size...
Interfaces and types that MONC components must specify.
logical function is_field_mse(name)
Determines whether a specific published field is a mse field.
real(kind=default_precision), dimension(:), allocatable wu_mse
real(kind=default_precision), dimension(:), allocatable th_buoyancy
Defined the local grid, i.e. the grid held on this process after decomposition.
real(kind=default_precision), dimension(:), allocatable ws_qt
subroutine populate_field_names()
Populates the published field names in the appropriate map.
Collection data structures.
real(kind=default_precision), dimension(:), allocatable sqt_qt
subroutine finalisation_callback(current_state)
Finalisation call back.
subroutine compute_uw_vw_for_column(current_state)
Computes the uw uv diagnostics for a specific column.
real(kind=default_precision), dimension(:,:), allocatable q_flux_values
real(kind=default_precision), dimension(:), allocatable w_thetal_advection
real(kind=default_precision), dimension(:,:), allocatable q_buoyancy
real(kind=default_precision) mflux
subroutine timestep_callback(current_state)
Timestep call back, this will deduce the diagnostics for the current (non halo) column.
real(kind=default_precision), dimension(:), allocatable thetal_thetal
real(kind=default_precision), dimension(:), allocatable wthetal_thetal
real(kind=default_precision), dimension(:,:), allocatable q_diff
real(kind=default_precision), dimension(:), allocatable u_mse_advection
logical some_prognostic_budget_diagnostics_enabled
real(kind=default_precision), dimension(:), allocatable w_mse_advection
logical function is_field_prognostic_budget(name)
Determines whether a specific published field is a uu, vv or ww field.
real(kind=default_precision), dimension(:), allocatable wu_u
real(kind=default_precision), dimension(:), allocatable ww_advection
real(kind=default_precision) function, public options_get_real(options_database, key, index)
Retrieves a real value from the database that matches the provided key.
real(kind=default_precision), dimension(:), allocatable v_mse
Functionality to support the different types of grid and abstraction between global grids and local o...
real(kind=default_precision), dimension(:), allocatable uw_viscosity
logical function is_field_q_flux(name)
Determines whether a specific published field is a q flux field.
integer function, public options_get_integer(options_database, key, index)
Retrieves an integer value from the database that matches the provided key.
subroutine compute_mse_for_column(current_state)
Computes the mse diagnostics for a specific column. For now we are assuming mse is the same as theta...
logical function is_field_uw_vw(name)
Determines whether a specific published field is a uw or uv field.
Manages the options database. Contains administration functions and deduce runtime options from the c...
real(kind=default_precision), dimension(:), allocatable ww_thetal
real(kind=default_precision), dimension(:), allocatable qt_qt
subroutine initialise_prognostic_budget_diagnostics(current_state)
Initialises the prognostic (uu, vv, ww) budget diagnostics.
real(kind=default_precision), dimension(:), allocatable v_thetal
real(kind=default_precision), dimension(:), allocatable mse_mse_diffusion
real(kind=default_precision), dimension(:,:), allocatable q_gradient
subroutine clear_q_fluxes()
Clears the Q flux diagnostics, called at the start of a timestep.
type(hashmap_type) scalar_fields
real(kind=default_precision), dimension(:), allocatable wqt_qt
subroutine field_value_retrieval_callback(current_state, name, field_value)
Field value retrieval callback, this returns the value of a specific published field.
subroutine compute_theta_flux_for_column(current_state)
Computes the heat flux diagnostics for a specific column.
real(kind=default_precision), dimension(:), allocatable v_qt_advection
real(kind=default_precision), dimension(:), allocatable uw_tendency
real(kind=default_precision), dimension(:), allocatable th_tendency
real(kind=default_precision), dimension(:), allocatable w_qt_buoyancy
Flux budget component which produces diagnostic data for the flux aspects of the model.
Determines whether or not a map contains a specific key.
logical some_q_flux_diagnostics_enabled
real(kind=default_precision), dimension(:), allocatable vs_qt
real(kind=default_precision), dimension(:), allocatable uw_buoyancy
subroutine clear_mse()
Clears the mse diagnostics.
real(kind=default_precision), dimension(:), allocatable thetal_thetal_diffusion
real(kind=default_precision), dimension(:), allocatable wv_thetal
real(kind=default_precision), dimension(:), allocatable qt_qt_diffusion
real(kind=default_precision), dimension(:), allocatable u_thetal_advection
real(kind=default_precision), dimension(:), allocatable wmse_mse
real(kind=default_precision), dimension(:), allocatable vw_tendency
subroutine compute_prognostic_budgets_for_column(current_state)
Computes the prognostic (uu, vv, ww) budgets for a specific column.
The model state which represents the current state of a run.
real(kind=default_precision), dimension(:), allocatable uu_advection
integer, parameter, public y_index
real(kind=default_precision), dimension(:), allocatable vw_w
real(kind=default_precision), dimension(:), allocatable wv_v
real(kind=default_precision), dimension(:), allocatable tv_sv
type(hashmap_type) uw_vw_fields
integer diagnostic_generation_frequency
real(kind=default_precision), dimension(:), allocatable sthetal_thetal
real(kind=default_precision), dimension(:), allocatable ww_viscosity
type(hashmap_type) heat_flux_fields
real(kind=default_precision), dimension(:), allocatable tw_sw
logical function is_field_heat_flux(name)
Determines whether a specific published field is a heat flux field.
real(kind=default_precision), dimension(:), allocatable vw_buoyancy
integer, parameter, public x_index
subroutine compute_qt_for_column(current_state)
Computes the qt diagnostics for a specific column. For now we are assuming qt is the same as theta...
real(kind=default_precision), dimension(:), allocatable wv_mse
real(kind=default_precision), dimension(:), allocatable vs_mse
logical some_qt_diagnostics_enabled
logical function is_field_thetal(name)
Determines whether a specific published field is a thetal field.
logical function get_published_field_enabled_state(collection, field_name)
Retrieves whether a published field is enabled or not.
integer, parameter, public component_double_data_type