This subroutine calculates and writes the total model energy, the energy and mass of each layer, and other globally integrated physical quantities.
298 type(ocean_grid_type),
intent(in) :: g
299 type(verticalgrid_type),
intent(in) :: gv
300 type(unit_scale_type),
intent(in) :: us
301 real,
dimension(SZIB_(G),SZJ_(G),SZK_(G)), &
303 real,
dimension(SZI_(G),SZJB_(G),SZK_(G)), &
305 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)), &
307 type(thermo_var_ptrs),
intent(in) :: tv
309 type(time_type),
intent(in) :: day
310 integer,
intent(in) :: n
312 type(sum_output_cs),
pointer :: cs
314 type(tracer_flow_control_cs), &
315 optional,
pointer :: tracer_csp
316 type(ocean_obc_type), &
317 optional,
pointer :: obc
318 type(time_type),
optional,
intent(in) :: dt_forcing
320 real :: eta(szi_(g),szj_(g),szk_(g)+1)
321 real :: areatm(szi_(g),szj_(g))
323 real :: pe(szk_(g)+1)
326 real :: z_0ape(szk_(g)+1)
328 real :: h_0ape(szk_(g)+1)
333 real :: vol_lay(szk_(g))
335 real :: mass_lay(szk_(g))
379 type(
efp_type),
dimension(5) :: efp_list
385 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)) :: &
387 real,
dimension(SZI_(G),SZJ_(G),SZK_(G)+1) :: &
389 real,
dimension(SZI_(G),SZJ_(G)) :: &
392 real :: ke_scale_factor
394 real :: pe_scale_factor
396 integer :: num_nc_fields
398 integer :: i, j, k, is, ie, js, je, ns, nz, m, isq, ieq, jsq, jeq, isr, ier, jsr, jer
399 integer :: l, lbelow, labove
402 integer :: start_of_day, num_days
404 character(len=240) :: energypath_nc
405 character(len=200) :: mesg
406 character(len=32) :: mesg_intro, time_units, day_str, n_str, date_str
407 logical :: date_stamped
408 type(time_type) :: dt_force
409 real :: tr_stocks(max_fields_)
410 real :: tr_min(max_fields_), tr_max(max_fields_)
411 real :: tr_min_x(max_fields_), tr_min_y(max_fields_), tr_min_z(max_fields_)
412 real :: tr_max_x(max_fields_), tr_max_y(max_fields_), tr_max_z(max_fields_)
413 logical :: tr_minmax_got(max_fields_) = .false.
414 character(len=40),
dimension(MAX_FIELDS_) :: &
416 integer :: ntr_stocks
417 integer :: iyear, imonth, iday, ihour, iminute, isecond, itick
418 logical :: local_open_bc
419 type(obc_segment_type),
pointer :: segment => null()
422 type(vardesc) :: vars(num_fields+max_fields_)
425 dt_force = set_time(seconds=2) ;
if (
present(dt_forcing)) dt_force = dt_forcing
426 if (cs%previous_calls == 0)
then 427 if (cs%energysave_geometric)
then 428 if (cs%energysavedays_geometric < cs%energysavedays)
then 429 cs%write_energy_time = day + cs%energysavedays_geometric
430 cs%geometric_end_time = cs%Start_time + cs%energysavedays * &
431 (1 + (day - cs%Start_time) / cs%energysavedays)
433 cs%write_energy_time = cs%Start_time + cs%energysavedays * &
434 (1 + (day - cs%Start_time) / cs%energysavedays)
437 cs%write_energy_time = cs%Start_time + cs%energysavedays * &
438 (1 + (day - cs%Start_time) / cs%energysavedays)
440 elseif (day + (dt_force/2) <= cs%write_energy_time)
then 443 if (cs%energysave_geometric)
then 444 if (cs%write_energy_time + cs%energysavedays_geometric >= &
445 cs%geometric_end_time)
then 446 cs%write_energy_time = cs%geometric_end_time
447 cs%energysave_geometric = .false.
449 cs%write_energy_time = cs%write_energy_time + cs%energysavedays_geometric
451 cs%energysavedays_geometric = cs%energysavedays_geometric*2
453 cs%write_energy_time = cs%write_energy_time + cs%energysavedays
458 if (.not.cs%use_temperature) num_nc_fields = 11
459 vars(1) = var_desc(
"Ntrunc",
"Nondim",
"Number of Velocity Truncations",
'1',
'1')
460 vars(2) = var_desc(
"En",
"Joules",
"Total Energy",
'1',
'1')
461 vars(3) = var_desc(
"APE",
"Joules",
"Total Interface APE",
'1',
'i')
462 vars(4) = var_desc(
"KE",
"Joules",
"Total Layer KE",
'1',
'L')
463 vars(5) = var_desc(
"H0",
"meter",
"Zero APE Depth of Interface",
'1',
'i')
464 vars(6) = var_desc(
"Mass_lay",
"kg",
"Total Layer Mass",
'1',
'L')
465 vars(7) = var_desc(
"Mass",
"kg",
"Total Mass",
'1',
'1')
466 vars(8) = var_desc(
"Mass_chg",
"kg",
"Total Mass Change between Entries",
'1',
'1')
467 vars(9) = var_desc(
"Mass_anom",
"kg",
"Anomalous Total Mass Change",
'1',
'1')
468 vars(10) = var_desc(
"max_CFL_trans",
"Nondim",
"Maximum finite-volume CFL",
'1',
'1')
469 vars(11) = var_desc(
"max_CFL_lin",
"Nondim",
"Maximum finite-difference CFL",
'1',
'1')
470 if (cs%use_temperature)
then 471 vars(12) = var_desc(
"Salt",
"kg",
"Total Salt",
'1',
'1')
472 vars(13) = var_desc(
"Salt_chg",
"kg",
"Total Salt Change between Entries",
'1',
'1')
473 vars(14) = var_desc(
"Salt_anom",
"kg",
"Anomalous Total Salt Change",
'1',
'1')
474 vars(15) = var_desc(
"Heat",
"Joules",
"Total Heat",
'1',
'1')
475 vars(16) = var_desc(
"Heat_chg",
"Joules",
"Total Heat Change between Entries",
'1',
'1')
476 vars(17) = var_desc(
"Heat_anom",
"Joules",
"Anomalous Total Heat Change",
'1',
'1')
479 local_open_bc = .false.
480 if (
present(obc))
then ;
if (
associated(obc))
then 481 local_open_bc = (obc%open_u_BCs_exist_globally .or. obc%open_v_BCs_exist_globally)
484 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec ; nz = g%ke
485 isq = g%IscB ; ieq = g%IecB ; jsq = g%JscB ; jeq = g%JecB
486 isr = is - (g%isd-1) ; ier = ie - (g%isd-1) ; jsr = js - (g%jsd-1) ; jer = je - (g%jsd-1)
489 hl2_to_kg = gv%H_to_kg_m2*us%L_to_m**2
491 if (.not.
associated(cs))
call mom_error(fatal, &
492 "write_energy: Module must be initialized before it is used.")
494 do j=js,je ;
do i=is,ie
495 areatm(i,j) = g%mask2dT(i,j)*g%areaT(i,j)
498 if (gv%Boussinesq)
then 500 do k=1,nz ;
do j=js,je ;
do i=is,ie
501 tmp1(i,j,k) = h(i,j,k) * (hl2_to_kg*areatm(i,j))
502 enddo ;
enddo ;
enddo 507 if (local_open_bc)
then 508 do ns=1, obc%number_of_segments
509 segment => obc%segment(ns)
510 if (.not. segment%on_pe .or. segment%specified) cycle
511 i=segment%HI%IsdB ; j=segment%HI%JsdB
512 if (segment%direction == obc_direction_e)
then 513 do k=1,nz ;
do j=segment%HI%jsd,segment%HI%jed
516 elseif (segment%direction == obc_direction_w)
then 517 do k=1,nz ;
do j=segment%HI%jsd,segment%HI%jed
520 elseif (segment%direction == obc_direction_n)
then 521 do k=1,nz ;
do i=segment%HI%isd,segment%HI%ied
524 elseif (segment%direction == obc_direction_s)
then 525 do k=1,nz ;
do i=segment%HI%isd,segment%HI%ied
532 mass_tot =
reproducing_sum(tmp1, isr, ier, jsr, jer, sums=mass_lay, efp_sum=mass_efp)
533 do k=1,nz ; vol_lay(k) = (us%m_to_L**2*gv%H_to_Z/gv%H_to_kg_m2)*mass_lay(k) ;
enddo 536 if (cs%do_APE_calc)
then 537 do k=1,nz ;
do j=js,je ;
do i=is,ie
538 tmp1(i,j,k) = hl2_to_kg * h(i,j,k) * areatm(i,j)
539 enddo ;
enddo ;
enddo 540 mass_tot =
reproducing_sum(tmp1, isr, ier, jsr, jer, sums=mass_lay, efp_sum=mass_efp)
542 call find_eta(h, tv, g, gv, us, eta)
543 do k=1,nz ;
do j=js,je ;
do i=is,ie
544 tmp1(i,j,k) = us%Z_to_m*us%L_to_m**2*(eta(i,j,k)-eta(i,j,k+1)) * areatm(i,j)
545 enddo ;
enddo ;
enddo 547 do k=1,nz ; vol_lay(k) = us%m_to_Z*us%m_to_L**2 * vol_lay(k) ;
enddo 549 do k=1,nz ;
do j=js,je ;
do i=is,ie
550 tmp1(i,j,k) = hl2_to_kg * h(i,j,k) * areatm(i,j)
551 enddo ;
enddo ;
enddo 552 mass_tot =
reproducing_sum(tmp1, isr, ier, jsr, jer, sums=mass_lay, efp_sum=mass_efp)
553 do k=1,nz ; vol_lay(k) = us%m_to_Z*us%m_to_L**2*us%kg_m3_to_R * (mass_lay(k) / gv%Rho0) ;
enddo 558 if (
present(tracer_csp))
then 559 call call_tracer_stocks(h, tr_stocks, g, gv, tracer_csp, stock_names=tr_names, &
560 stock_units=tr_units, num_stocks=ntr_stocks,&
561 got_min_max=tr_minmax_got, global_min=tr_min, global_max=tr_max, &
562 xgmin=tr_min_x, ygmin=tr_min_y, zgmin=tr_min_z,&
563 xgmax=tr_max_x, ygmax=tr_max_y, zgmax=tr_max_z)
564 if (ntr_stocks > 0)
then 566 vars(num_nc_fields+m) = var_desc(tr_names(m), units=tr_units(m), &
567 longname=tr_names(m), hor_grid=
'1', z_grid=
'1')
569 num_nc_fields = num_nc_fields + ntr_stocks
573 if (cs%previous_calls == 0)
then 575 cs%mass_prev_EFP = mass_efp
576 cs%fresh_water_in_EFP = real_to_efp(0.0)
577 if (cs%use_temperature)
then 578 cs%net_salt_in_EFP = real_to_efp(0.0) ; cs%net_heat_in_EFP = real_to_efp(0.0)
582 if (is_root_pe())
then 583 if (day > cs%Start_time)
then 584 call open_file(cs%fileenergy_ascii, trim(cs%energyfile), &
585 action=append_file, form=ascii_file, nohdrs=.true.)
587 call open_file(cs%fileenergy_ascii, trim(cs%energyfile), &
588 action=writeonly_file, form=ascii_file, nohdrs=.true.)
589 if (abs(cs%timeunit - 86400.0) < 1.0)
then 590 if (cs%use_temperature)
then 591 write(cs%fileenergy_ascii,
'(" Step,",7x,"Day, Truncs, & 592 &Energy/Mass, Maximum CFL, Mean Sea Level, Total Mass, Mean Salin, & 593 &Mean Temp, Frac Mass Err, Salin Err, Temp Err")')
594 write(cs%fileenergy_ascii,
'(12x,"[days]",17x,"[m2 s-2]",11x,"[Nondim]",7x,"[m]",13x,& 595 &"[kg]",9x,"[PSU]",6x,"[degC]",7x,"[Nondim]",8x,"[PSU]",8x,"[degC]")')
597 write(cs%fileenergy_ascii,
'(" Step,",7x,"Day, Truncs, & 598 &Energy/Mass, Maximum CFL, Mean sea level, Total Mass, Frac Mass Err")')
599 write(cs%fileenergy_ascii,
'(12x,"[days]",17x,"[m2 s-2]",11x,"[Nondim]",8x,"[m]",13x,& 600 &"[kg]",11x,"[Nondim]")')
603 if ((cs%timeunit >= 0.99) .and. (cs%timeunit < 1.01))
then 604 time_units =
" [seconds] " 605 elseif ((cs%timeunit >= 3599.0) .and. (cs%timeunit < 3601.0))
then 606 time_units =
" [hours] " 607 elseif ((cs%timeunit >= 86399.0) .and. (cs%timeunit < 86401.0))
then 608 time_units =
" [days] " 609 elseif ((cs%timeunit >= 3.0e7) .and. (cs%timeunit < 3.2e7))
then 610 time_units =
" [years] " 612 write(time_units,
'(9x,"[",es8.2," s] ")') cs%timeunit
615 if (cs%use_temperature)
then 616 write(cs%fileenergy_ascii,
'(" Step,",7x,"Time, Truncs, & 617 &Energy/Mass, Maximum CFL, Mean Sea Level, Total Mass, Mean Salin, & 618 &Mean Temp, Frac Mass Err, Salin Err, Temp Err")')
619 write(cs%fileenergy_ascii,
'(A25,10x,"[m2 s-2]",11x,"[Nondim]",7x,"[m]",13x,& 620 &"[kg]",9x,"[PSU]",6x,"[degC]",7x,"[Nondim]",8x,"[PSU]",6x,& 621 &"[degC]")') time_units
623 write(cs%fileenergy_ascii,
'(" Step,",7x,"Time, Truncs, & 624 &Energy/Mass, Maximum CFL, Mean sea level, Total Mass, Frac Mass Err")')
625 write(cs%fileenergy_ascii,
'(A25,10x,"[m2 s-2]",11x,"[Nondim]",8x,"[m]",13x,& 626 &"[kg]",11x,"[Nondim]")') time_units
632 energypath_nc = trim(cs%energyfile) //
".nc" 633 if (day > cs%Start_time)
then 634 call reopen_file(cs%fileenergy_nc, trim(energypath_nc), vars, &
635 num_nc_fields, cs%fields, single_file, cs%timeunit, &
638 call create_file(cs%fileenergy_nc, trim(energypath_nc), vars, &
639 num_nc_fields, cs%fields, single_file, cs%timeunit, &
644 if (cs%do_APE_calc)
then 645 lbelow = 1 ; volbelow = 0.0
647 volbelow = volbelow + vol_lay(k)
648 if ((volbelow >= cs%DL(cs%lH(k))%vol_below) .and. &
649 (volbelow < cs%DL(cs%lH(k)+1)%vol_below))
then 652 labove=cs%list_size+1
653 l = (labove + lbelow) / 2
654 do while (l > lbelow)
655 if (volbelow < cs%DL(l)%vol_below)
then ; labove = l
656 else ; lbelow = l ;
endif 657 l = (labove + lbelow) / 2
662 z_0ape(k) = cs%DL(l)%depth - (volbelow - cs%DL(l)%vol_below) / cs%DL(l)%area
664 z_0ape(nz+1) = cs%DL(2)%depth
669 pe_scale_factor = us%RZ_to_kg_m2*us%L_to_m**2*us%L_T_to_m_s**2
671 if (gv%Boussinesq)
then 672 do j=js,je ;
do i=is,ie
675 hbelow = hbelow + h(i,j,k) * gv%H_to_Z
676 hint = z_0ape(k) + (hbelow - g%bathyT(i,j))
677 hbot = z_0ape(k) - g%bathyT(i,j)
678 hbot = (hbot + abs(hbot)) * 0.5
679 pe_pt(i,j,k) = (0.5 * pe_scale_factor * areatm(i,j)) * (gv%Rho0*gv%g_prime(k)) * &
680 (hint * hint - hbot * hbot)
684 do j=js,je ;
do i=is,ie
686 hint = z_0ape(k) + eta(i,j,k)
687 hbot = max(z_0ape(k) - g%bathyT(i,j), 0.0)
688 pe_pt(i,j,k) = (0.5 * pe_scale_factor * areatm(i,j) * (gv%Rho0*gv%g_prime(k))) * &
689 (hint * hint - hbot * hbot)
695 do k=1,nz+1 ; h_0ape(k) = us%Z_to_m*z_0ape(k) ;
enddo 698 do k=1,nz+1 ; pe(k) = 0.0 ; h_0ape(k) = 0.0 ;
enddo 702 ke_scale_factor = hl2_to_kg*us%L_T_to_m_s**2
704 do k=1,nz ;
do j=js,je ;
do i=is,ie
705 tmp1(i,j,k) = (0.25 * ke_scale_factor * (areatm(i,j) * h(i,j,k))) * &
706 ((u(i-1,j,k)**2 + u(i,j,k)**2) + (v(i,j-1,k)**2 + v(i,j,k)**2))
707 enddo ;
enddo ;
enddo 711 toten = ke_tot + pe_tot
713 salt = 0.0 ; heat = 0.0
714 if (cs%use_temperature)
then 715 temp_int(:,:) = 0.0 ; salt_int(:,:) = 0.0
716 do k=1,nz ;
do j=js,je ;
do i=is,ie
717 salt_int(i,j) = salt_int(i,j) + tv%S(i,j,k) * &
718 (h(i,j,k)*(hl2_to_kg * areatm(i,j)))
719 temp_int(i,j) = temp_int(i,j) + (us%Q_to_J_kg*tv%C_p * tv%T(i,j,k)) * &
720 (h(i,j,k)*(hl2_to_kg * areatm(i,j)))
721 enddo ;
enddo ;
enddo 726 efp_list(1) = salt_efp ; efp_list(2) = heat_efp ; efp_list(3) = cs%fresh_water_in_EFP
727 efp_list(4) = cs%net_salt_in_EFP ; efp_list(5) = cs%net_heat_in_EFP
730 salt_efp = efp_list(1) ; heat_efp = efp_list(2) ; cs%fresh_water_in_EFP = efp_list(3)
731 cs%net_salt_in_EFP = efp_list(4) ; cs%net_heat_in_EFP = efp_list(5)
733 salt = efp_to_real(salt_efp)
734 heat = efp_to_real(heat_efp)
741 do k=1,nz ;
do j=js,je ;
do i=isq,ieq
742 cfl_iarea = g%IareaT(i,j)
743 if (u(i,j,k) < 0.0) &
744 cfl_iarea = g%IareaT(i+1,j)
746 cfl_trans = abs(u(i,j,k) * cs%dt_in_T) * (g%dy_Cu(i,j) * cfl_iarea)
747 cfl_lin = abs(u(i,j,k) * cs%dt_in_T) * g%IdxCu(i,j)
748 max_cfl(1) = max(max_cfl(1), cfl_trans)
749 max_cfl(2) = max(max_cfl(2), cfl_lin)
750 enddo ;
enddo ;
enddo 751 do k=1,nz ;
do j=jsq,jeq ;
do i=is,ie
752 cfl_iarea = g%IareaT(i,j)
753 if (v(i,j,k) < 0.0) &
754 cfl_iarea = g%IareaT(i,j+1)
756 cfl_trans = abs(v(i,j,k) * cs%dt_in_T) * (g%dx_Cv(i,j) * cfl_iarea)
757 cfl_lin = abs(v(i,j,k) * cs%dt_in_T) * g%IdyCv(i,j)
758 max_cfl(1) = max(max_cfl(1), cfl_trans)
759 max_cfl(2) = max(max_cfl(2), cfl_lin)
760 enddo ;
enddo ;
enddo 762 call sum_across_pes(cs%ntrunc)
766 if (ntr_stocks > 0)
call sum_across_pes(tr_stocks,ntr_stocks)
768 call max_across_pes(max_cfl, 2)
769 irho0 = 1.0 / (us%R_to_kg_m3*gv%Rho0)
771 if (cs%use_temperature)
then 772 if (cs%previous_calls == 0)
then 773 cs%salt_prev_EFP = salt_efp ; cs%heat_prev_EFP = heat_efp
775 salt_chg_efp = salt_efp - cs%salt_prev_EFP
776 salt_chg = efp_to_real(salt_chg_efp)
777 salt_anom_efp = salt_chg_efp - cs%net_salt_in_EFP
778 salt_anom = efp_to_real(salt_anom_efp)
779 heat_chg_efp = heat_efp - cs%heat_prev_EFP
780 heat_chg = efp_to_real(heat_chg_efp)
781 heat_anom_efp = heat_chg_efp - cs%net_heat_in_EFP
782 heat_anom = efp_to_real(heat_anom_efp)
785 mass_chg_efp = mass_efp - cs%mass_prev_EFP
786 mass_anom_efp = mass_chg_efp - cs%fresh_water_in_EFP
787 mass_anom = efp_to_real(mass_anom_efp)
788 if (cs%use_temperature .and. .not.gv%Boussinesq)
then 790 mass_anom = mass_anom - 0.001*efp_to_real(cs%net_salt_in_EFP)
792 mass_chg = efp_to_real(mass_chg_efp)
794 if (cs%use_temperature)
then 795 salin = salt / mass_tot ; salin_anom = salt_anom / mass_tot
797 temp = heat / (mass_tot*us%Q_to_J_kg*tv%C_p) ; temp_anom = heat_anom / (mass_tot*us%Q_to_J_kg*tv%C_p)
799 en_mass = toten / mass_tot
801 call get_time(day, start_of_day, num_days)
802 date_stamped = (cs%date_stamped_output .and. (get_calendar_type() /= no_calendar))
804 call get_date(day, iyear, imonth, iday, ihour, iminute, isecond, itick)
805 if (abs(cs%timeunit - 86400.0) < 1.0)
then 806 reday =
REAL(num_days)+ (
REAL(start_of_day)/86400.0)
807 mesg_intro =
"MOM Day " 809 reday =
REAL(num_days)*(86400.0/cs%timeunit) + &
810 REAL(start_of_day)/abs(cs%timeunit)
811 mesg_intro =
"MOM Time " 813 if (reday < 1.0e8)
then ;
write(day_str,
'(F12.3)') reday
814 elseif (reday < 1.0e11)
then ;
write(day_str,
'(F15.3)') reday
815 else ;
write(day_str,
'(ES15.9)') reday ;
endif 817 if (n < 1000000)
then ;
write(n_str,
'(I6)') n
818 elseif (n < 10000000)
then ;
write(n_str,
'(I7)') n
819 elseif (n < 100000000)
then ;
write(n_str,
'(I8)') n
820 else ;
write(n_str,
'(I10)') n ;
endif 822 if (date_stamped)
then 823 write(date_str,
'("MOM Date",i7,2("/",i2.2)," ",i2.2,2(":",i2.2))') &
824 iyear, imonth, iday, ihour, iminute, isecond
826 date_str = trim(mesg_intro)//trim(day_str)
829 if (is_root_pe())
then 830 if (cs%use_temperature)
then 831 write(*,
'(A," ",A,": En ",ES12.6, ", MaxCFL ", F8.5, ", Mass ", & 832 & ES18.12, ", Salt ", F15.11,", Temp ", F15.11)') &
833 trim(date_str), trim(n_str), en_mass, max_cfl(1), mass_tot, salin, temp
835 write(*,
'(A," ",A,": En ",ES12.6, ", MaxCFL ", F8.5, ", Mass ", & 837 trim(date_str), trim(n_str), en_mass, max_cfl(1), mass_tot
840 if (cs%use_temperature)
then 841 write(cs%fileenergy_ascii,
'(A,",",A,",", I6,", En ",ES22.16, & 842 &", CFL ", F8.5, ", SL ",& 843 &es11.4,", M ",ES11.5,", S",f8.4,", T",f8.4,& 844 &", Me ",ES9.2,", Se ",ES9.2,", Te ",ES9.2)') &
845 trim(n_str), trim(day_str), cs%ntrunc, en_mass, max_cfl(1), &
846 -h_0ape(1), mass_tot, salin, temp, mass_anom/mass_tot, salin_anom, &
849 write(cs%fileenergy_ascii,
'(A,",",A,",", I6,", En ",ES22.16, & 850 &", CFL ", F8.5, ", SL ",& 851 &ES11.4,", Mass ",ES11.5,", Me ",ES9.2)') &
852 trim(n_str), trim(day_str), cs%ntrunc, en_mass, max_cfl(1), &
853 -h_0ape(1), mass_tot, mass_anom/mass_tot
856 if (cs%ntrunc > 0)
then 857 write(*,
'(A," Energy/Mass:",ES12.5," Truncations ",I0)') &
858 trim(date_str), en_mass, cs%ntrunc
861 if (cs%write_stocks)
then 862 write(*,
'(" Total Energy: ",Z16.16,ES24.16)') toten, toten
863 write(*,
'(" Total Mass: ",ES24.16,", Change: ",ES24.16," Error: ",ES12.5," (",ES8.1,")")') &
864 mass_tot, mass_chg, mass_anom, mass_anom/mass_tot
865 if (cs%use_temperature)
then 867 write(*,
'(" Total Salt: ",ES24.16,", Change: ",ES24.16," Error: ",ES12.5)') &
868 salt*0.001, salt_chg*0.001, salt_anom*0.001
870 write(*,
'(" Total Salt: ",ES24.16,", Change: ",ES24.16," Error: ",ES12.5," (",ES8.1,")")') &
871 salt*0.001, salt_chg*0.001, salt_anom*0.001, salt_anom/salt
874 write(*,
'(" Total Heat: ",ES24.16,", Change: ",ES24.16," Error: ",ES12.5)') &
875 heat, heat_chg, heat_anom
877 write(*,
'(" Total Heat: ",ES24.16,", Change: ",ES24.16," Error: ",ES12.5," (",ES8.1,")")') &
878 heat, heat_chg, heat_anom, heat_anom/heat
883 write(*,
'(" Total ",a,": ",ES24.16,X,a)') &
884 trim(tr_names(m)), tr_stocks(m), trim(tr_units(m))
886 if (tr_minmax_got(m))
then 887 write(*,
'(64X,"Global Min:",ES24.16,X,"at: (", f7.2,","f7.2,","f8.2,")" )') &
888 tr_min(m),tr_min_x(m),tr_min_y(m),tr_min_z(m)
889 write(*,
'(64X,"Global Max:",ES24.16,X,"at: (", f7.2,","f7.2,","f8.2,")" )') &
890 tr_max(m),tr_max_x(m),tr_max_y(m),tr_max_z(m)
897 var =
real(cs%ntrunc)
898 call write_field(cs%fileenergy_nc, cs%fields(1), var, reday)
899 call write_field(cs%fileenergy_nc, cs%fields(2), toten, reday)
900 call write_field(cs%fileenergy_nc, cs%fields(3), pe, reday)
901 call write_field(cs%fileenergy_nc, cs%fields(4), ke, reday)
902 call write_field(cs%fileenergy_nc, cs%fields(5), h_0ape, reday)
903 call write_field(cs%fileenergy_nc, cs%fields(6), mass_lay, reday)
905 call write_field(cs%fileenergy_nc, cs%fields(7), mass_tot, reday)
906 call write_field(cs%fileenergy_nc, cs%fields(8), mass_chg, reday)
907 call write_field(cs%fileenergy_nc, cs%fields(9), mass_anom, reday)
908 call write_field(cs%fileenergy_nc, cs%fields(10), max_cfl(1), reday)
909 call write_field(cs%fileenergy_nc, cs%fields(11), max_cfl(2), reday)
910 if (cs%use_temperature)
then 911 call write_field(cs%fileenergy_nc, cs%fields(12), 0.001*salt, reday)
912 call write_field(cs%fileenergy_nc, cs%fields(13), 0.001*salt_chg, reday)
913 call write_field(cs%fileenergy_nc, cs%fields(14), 0.001*salt_anom, reday)
914 call write_field(cs%fileenergy_nc, cs%fields(15), heat, reday)
915 call write_field(cs%fileenergy_nc, cs%fields(16), heat_chg, reday)
916 call write_field(cs%fileenergy_nc, cs%fields(17), heat_anom, reday)
918 call write_field(cs%fileenergy_nc, cs%fields(17+m), tr_stocks(m), reday)
922 call write_field(cs%fileenergy_nc, cs%fields(11+m), tr_stocks(m), reday)
926 call flush_file(cs%fileenergy_nc)
929 if ((en_mass>cs%max_Energy) .or. &
930 ((en_mass>cs%max_Energy) .and. (en_mass<cs%max_Energy)))
then 931 write(mesg,
'("Energy per unit mass of ",ES11.4," exceeds ",ES11.4)') &
932 en_mass, cs%max_Energy
933 call mom_error(fatal, &
934 "write_energy : Excessive energy per unit mass or NaNs forced model termination.")
936 if (cs%ntrunc>cs%maxtrunc)
then 937 call mom_error(fatal,
"write_energy : Ocean velocity has been truncated too many times.")
940 cs%previous_calls = cs%previous_calls + 1
942 cs%mass_prev_EFP = mass_efp ; cs%fresh_water_in_EFP = real_to_efp(0.0)
943 if (cs%use_temperature)
then 944 cs%salt_prev_EFP = salt_efp ; cs%net_salt_in_EFP = real_to_efp(0.0)
945 cs%heat_prev_EFP = heat_efp ; cs%net_heat_in_EFP = real_to_efp(0.0)
Sum a value or 1-d array of values across processors, returning the sums in place.
Find an accurate and order-invariant sum of a distributed 2d or 3d field.
Find an accurate and order-invariant sum of a distributed 2d field, returning the result in the form ...
The Extended Fixed Point (EFP) type provides a public interface for doing sums and taking differences...