6 use iso_fortran_env
, only : int64
7 use mom_coms, only : sum_across_pes, pe_here, root_pe, num_pes, max_across_pes
15 use mom_io, only : create_file, fieldtype, flush_file, open_file, reopen_file
17 use mom_io, only : append_file, ascii_file, single_file, writeonly_file
19 use mom_open_boundary, only : obc_direction_e, obc_direction_w, obc_direction_n, obc_direction_s
20 use mom_time_manager, only : time_type, get_time, get_date, set_time,
operator(>)
21 use mom_time_manager, only :
operator(+),
operator(-),
operator(*),
operator(/)
22 use mom_time_manager, only :
operator(/=),
operator(<=),
operator(>=),
operator(<)
23 use mom_time_manager, only : get_calendar_type, time_type_to_real, no_calendar
28 use mpp_mod
, only : mpp_chksum
32 implicit none ;
private 34 #include <MOM_memory.h> 36 public write_energy, accumulate_net_input, mom_sum_output_init
43 integer,
parameter :: num_fields = 17
44 character (*),
parameter :: depth_chksum_attr =
"bathyT_checksum" 47 character (*),
parameter :: area_chksum_attr =
"mask2dT_areaT_checksum" 65 integer,
allocatable,
dimension(:) :: lh
68 logical :: do_ape_calc
71 logical :: read_depth_list
73 character(len=200) :: depth_list_file
74 real :: d_list_min_inc
76 logical :: require_depth_list_chksum
79 logical :: update_depth_list_chksum
82 logical :: use_temperature
97 type(time_type) :: energysavedays
99 type(time_type) :: energysavedays_geometric
103 logical :: energysave_geometric
105 type(time_type) :: write_energy_time
106 type(time_type) :: geometric_end_time
111 logical :: date_stamped_output
112 type(time_type) :: start_time
114 integer,
pointer :: ntrunc => null()
120 logical :: write_stocks
122 integer :: previous_calls = 0
123 integer :: prev_n = 0
124 integer :: fileenergy_nc
125 integer :: fileenergy_ascii
126 type(fieldtype),
dimension(NUM_FIELDS+MAX_FIELDS_) :: &
128 character(len=200) :: energyfile
134 subroutine mom_sum_output_init(G, US, param_file, directory, ntrnc, &
135 Input_start_time, CS)
140 character(len=*),
intent(in) :: directory
141 integer,
target,
intent(inout) :: ntrnc
144 type(time_type),
intent(in) :: Input_start_time
152 #include "version_variable.h" 153 character(len=40) :: mdl =
"MOM_sum_output" 154 character(len=200) :: energyfile
155 character(len=32) :: filename_appendix =
'' 157 if (
associated(cs))
then 158 call mom_error(warning,
"MOM_sum_output_init called with associated control structure.")
165 call get_param(param_file, mdl,
"CALCULATE_APE", cs%do_APE_calc, &
166 "If true, calculate the available potential energy of "//&
167 "the interfaces. Setting this to false reduces the "//&
168 "memory footprint of high-PE-count models dramatically.", &
170 call get_param(param_file, mdl,
"WRITE_STOCKS", cs%write_stocks, &
171 "If true, write the integrated tracer amounts to stdout "//&
172 "when the energy files are written.", default=.true.)
173 call get_param(param_file, mdl,
"ENABLE_THERMODYNAMICS", cs%use_temperature, &
174 "If true, Temperature and salinity are used as state "//&
175 "variables.", default=.true.)
176 call get_param(param_file, mdl,
"DT", cs%dt_in_T, &
177 "The (baroclinic) dynamics time step.", &
178 units=
"s", scale=us%s_to_T, fail_if_missing=.true.)
179 call get_param(param_file, mdl,
"MAXTRUNC", cs%maxtrunc, &
180 "The run will be stopped, and the day set to a very "//&
181 "large value if the velocity is truncated more than "//&
182 "MAXTRUNC times between energy saves. Set MAXTRUNC to 0 "//&
183 "to stop if there is any truncation of velocities.", &
184 units=
"truncations save_interval-1", default=0)
186 call get_param(param_file, mdl,
"MAX_ENERGY", cs%max_Energy, &
187 "The maximum permitted average energy per unit mass; the "//&
188 "model will be stopped if there is more energy than "//&
189 "this. If zero or negative, this is set to 10*MAXVEL^2.", &
190 units=
"m2 s-2", default=0.0)
191 if (cs%max_Energy <= 0.0)
then 192 call get_param(param_file, mdl,
"MAXVEL", maxvel, &
193 "The maximum velocity allowed before the velocity "//&
194 "components are truncated.", units=
"m s-1", default=3.0e8)
195 cs%max_Energy = 10.0 * maxvel**2
196 call log_param(param_file, mdl,
"MAX_ENERGY as used", cs%max_Energy)
199 call get_param(param_file, mdl,
"ENERGYFILE", energyfile, &
200 "The file to use to write the energies and globally "//&
201 "summed diagnostics.", default=
"ocean.stats")
204 call get_filename_appendix(filename_appendix)
205 if (len_trim(filename_appendix) > 0)
then 206 energyfile = trim(energyfile) //
'.'//trim(filename_appendix)
209 cs%energyfile = trim(slasher(directory))//trim(energyfile)
210 call log_param(param_file, mdl,
"output_path/ENERGYFILE", cs%energyfile)
212 cs%energyfile = trim(cs%energyfile)//
"."//trim(adjustl(statslabel))
215 call get_param(param_file, mdl,
"DATE_STAMPED_STDOUT", cs%date_stamped_output, &
216 "If true, use dates (not times) in messages to stdout", &
218 call get_param(param_file, mdl,
"TIMEUNIT", cs%Timeunit, &
219 "The time unit in seconds a number of input fields", &
220 units=
"s", default=86400.0)
221 if (cs%Timeunit < 0.0) cs%Timeunit = 86400.0
225 if (cs%do_APE_calc)
then 226 call get_param(param_file, mdl,
"READ_DEPTH_LIST", cs%read_depth_list, &
227 "Read the depth list from a file if it exists or "//&
228 "create that file otherwise.", default=.false.)
229 call get_param(param_file, mdl,
"DEPTH_LIST_MIN_INC", cs%D_list_min_inc, &
230 "The minimum increment between the depths of the "//&
231 "entries in the depth-list file.", &
232 units=
"m", default=1.0e-10, scale=us%m_to_Z)
233 if (cs%read_depth_list)
then 234 call get_param(param_file, mdl,
"DEPTH_LIST_FILE", cs%depth_list_file, &
235 "The name of the depth list file.", default=
"Depth_list.nc")
236 if (scan(cs%depth_list_file,
'/') == 0) &
237 cs%depth_list_file = trim(slasher(directory))//trim(cs%depth_list_file)
239 call get_param(param_file, mdl,
"REQUIRE_DEPTH_LIST_CHECKSUMS", &
240 cs%require_depth_list_chksum, &
241 "Require that matching checksums be in Depth_list.nc "//&
242 "when reading the file.", default=.true.)
243 if (.not. cs%require_depth_list_chksum) &
244 call get_param(param_file, mdl,
"UPDATE_DEPTH_LIST_CHECKSUMS", &
245 cs%update_depth_list_chksum, &
246 "Automatically update the Depth_list.nc file if the "//&
247 "checksums are missing or do not match current values.", &
251 allocate(cs%lH(g%ke))
252 call depth_list_setup(g, us, cs)
257 call get_param(param_file, mdl,
"TIMEUNIT", time_unit, &
258 "The time unit for ENERGYSAVEDAYS.", &
259 units=
"s", default=86400.0)
260 call get_param(param_file, mdl,
"ENERGYSAVEDAYS",cs%energysavedays, &
261 "The interval in units of TIMEUNIT between saves of the "//&
262 "energies of the run and other globally summed diagnostics.",&
263 default=set_time(0,days=1), timeunit=time_unit)
264 call get_param(param_file, mdl,
"ENERGYSAVEDAYS_GEOMETRIC",cs%energysavedays_geometric, &
265 "The starting interval in units of TIMEUNIT for the first call "//&
266 "to save the energies of the run and other globally summed diagnostics. "//&
267 "The interval increases by a factor of 2. after each call to write_energy.",&
268 default=set_time(seconds=0), timeunit=time_unit)
270 if ((time_type_to_real(cs%energysavedays_geometric) > 0.) .and. &
271 (cs%energysavedays_geometric < cs%energysavedays))
then 272 cs%energysave_geometric = .true.
274 cs%energysave_geometric = .false.
277 cs%Start_time = input_start_time
280 end subroutine mom_sum_output_init
283 subroutine mom_sum_output_end(CS)
286 if (
associated(cs))
then 287 if (cs%do_APE_calc)
then 288 deallocate(cs%lH, cs%DL)
293 end subroutine mom_sum_output_end
297 subroutine write_energy(u, v, h, tv, day, n, G, GV, US, CS, tracer_CSp, OBC, dt_forcing)
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)), &
309 type(time_type),
intent(in) :: day
310 integer,
intent(in) :: n
315 optional,
pointer :: tracer_CSp
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
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)
948 end subroutine write_energy
952 subroutine accumulate_net_input(fluxes, sfc_state, tv, dt, G, US, CS)
955 type(
surface),
intent(in) :: sfc_state
959 real,
intent(in) :: dt
965 real,
dimension(SZI_(G),SZJ_(G)) :: &
989 integer :: i, j, is, ie, js, je, isr, ier, jsr, jer
991 is = g%isc ; ie = g%iec ; js = g%jsc ; je = g%jec
993 rzl2_to_kg = us%L_to_m**2*us%RZ_to_kg_m2
994 qrzl2_to_j = rzl2_to_kg*us%Q_to_J_kg
997 if (
associated(fluxes%evap))
then 998 if (
associated(fluxes%lprec) .and.
associated(fluxes%fprec))
then 999 do j=js,je ;
do i=is,ie
1000 fw_in(i,j) = rzl2_to_kg * dt*g%areaT(i,j)*(fluxes%evap(i,j) + &
1001 (((fluxes%lprec(i,j) + fluxes%vprec(i,j)) + fluxes%lrunoff(i,j)) + &
1002 (fluxes%fprec(i,j) + fluxes%frunoff(i,j))))
1005 call mom_error(warning, &
1006 "accumulate_net_input called with associated evap field, but no precip field.")
1010 if (
associated(fluxes%seaice_melt))
then ;
do j=js,je ;
do i=is,ie
1011 fw_in(i,j) = fw_in(i,j) + rzl2_to_kg*dt * &
1012 g%areaT(i,j) * fluxes%seaice_melt(i,j)
1013 enddo ;
enddo ;
endif 1015 salt_in(:,:) = 0.0 ; heat_in(:,:) = 0.0
1016 if (cs%use_temperature)
then 1018 if (
associated(fluxes%sw))
then ;
do j=js,je ;
do i=is,ie
1019 heat_in(i,j) = heat_in(i,j) + dt * qrzl2_to_j * g%areaT(i,j) * (fluxes%sw(i,j) + &
1020 (fluxes%lw(i,j) + (fluxes%latent(i,j) + fluxes%sens(i,j))))
1021 enddo ;
enddo ;
endif 1023 if (
associated(fluxes%seaice_melt_heat))
then ;
do j=js,je ;
do i=is,ie
1024 heat_in(i,j) = heat_in(i,j) + dt * qrzl2_to_j * g%areaT(i,j) * &
1025 fluxes%seaice_melt_heat(i,j)
1026 enddo ;
enddo ;
endif 1039 if (
associated(tv%TempxPmE))
then 1040 do j=js,je ;
do i=is,ie
1041 heat_in(i,j) = heat_in(i,j) + (fluxes%C_p * qrzl2_to_j*g%areaT(i,j)) * tv%TempxPmE(i,j)
1043 elseif (
associated(fluxes%evap))
then 1044 do j=js,je ;
do i=is,ie
1045 heat_in(i,j) = heat_in(i,j) + (us%Q_to_J_kg*fluxes%C_p * sfc_state%SST(i,j)) * fw_in(i,j)
1050 if (
associated(tv%internal_heat))
then 1051 do j=js,je ;
do i=is,ie
1052 heat_in(i,j) = heat_in(i,j) + (fluxes%C_p * qrzl2_to_j*g%areaT(i,j)) * tv%internal_heat(i,j)
1055 if (
associated(tv%frazil))
then ;
do j=js,je ;
do i=is,ie
1056 heat_in(i,j) = heat_in(i,j) + qrzl2_to_j * g%areaT(i,j) * tv%frazil(i,j)
1057 enddo ;
enddo ;
endif 1058 if (
associated(fluxes%heat_added))
then ;
do j=js,je ;
do i=is,ie
1059 heat_in(i,j) = heat_in(i,j) + qrzl2_to_j * dt*g%areaT(i,j) * fluxes%heat_added(i,j)
1060 enddo ;
enddo ;
endif 1065 if (
associated(fluxes%salt_flux))
then ;
do j=js,je ;
do i=is,ie
1067 salt_in(i,j) = rzl2_to_kg * dt * &
1068 g%areaT(i,j)*(1000.0*fluxes%salt_flux(i,j))
1069 enddo ;
enddo ;
endif 1072 if ((cs%use_temperature) .or.
associated(fluxes%lprec) .or. &
1073 associated(fluxes%evap))
then 1076 isr = is - (g%isd-1) ; ier = ie - (g%isd-1) ; jsr = js - (g%jsd-1) ; jer = je - (g%jsd-1)
1081 cs%fresh_water_in_EFP = cs%fresh_water_in_EFP + fw_in_efp
1082 cs%net_salt_in_EFP = cs%net_salt_in_EFP + salt_in_efp
1083 cs%net_heat_in_EFP = cs%net_heat_in_EFP + heat_in_efp
1086 end subroutine accumulate_net_input
1092 subroutine depth_list_setup(G, US, CS)
1100 if (cs%read_depth_list)
then 1102 call read_depth_list(g, us, cs, cs%depth_list_file)
1104 if (is_root_pe())
call mom_error(warning,
"depth_list_setup: "// &
1105 trim(cs%depth_list_file)//
" does not exist. Creating a new file.")
1106 call create_depth_list(g, cs)
1108 call write_depth_list(g, us, cs, cs%depth_list_file, cs%list_size+1)
1111 call create_depth_list(g, cs)
1115 cs%lH(k) = cs%list_size
1118 end subroutine depth_list_setup
1122 subroutine create_depth_list(G, CS)
1127 real,
dimension(G%Domain%niglobal*G%Domain%njglobal + 1) :: &
1128 Dlist, & !< The global list of bottom depths [Z ~> m].
1130 integer,
dimension(G%Domain%niglobal*G%Domain%njglobal+1) :: &
1137 logical :: add_to_list
1139 integer :: ir, indxt
1140 integer :: mls, list_size
1141 integer :: list_pos, i_global, j_global
1142 integer :: i, j, k, kl
1144 mls = g%Domain%niglobal*g%Domain%njglobal
1149 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1151 j_global = j + g%jdg_offset - (g%jsg-1)
1152 i_global = i + g%idg_offset - (g%isg-1)
1154 list_pos = (j_global-1)*g%Domain%niglobal + i_global
1155 dlist(list_pos) = g%bathyT(i,j)
1156 arealist(list_pos) = g%mask2dT(i,j) * g%areaT(i,j)
1160 call sum_across_pes(dlist, mls+1)
1161 call sum_across_pes(arealist, mls+1)
1163 do j=1,mls+1 ; indx2(j) = j ;
enddo 1164 k = mls / 2 + 1 ; ir = mls
1173 indx2(ir) = indx2(1)
1175 if (ir == 1)
then ; indx2(1) = indxt ;
exit ;
endif 1178 do ;
if (j > ir)
exit 1179 if (j < ir .AND. dlist(indx2(j)) < dlist(indx2(j+1))) j = j + 1
1180 if (dnow < dlist(indx2(j)))
then ; indx2(i) = indx2(j) ; i = j ; j = j + i
1181 else ; j = ir+1 ;
endif 1193 d_list_prev = dlist(indx2(mls))
1196 if (dlist(indx2(k)) < d_list_prev-cs%D_list_min_inc)
then 1197 list_size = list_size + 1
1198 d_list_prev = dlist(indx2(k))
1202 cs%list_size = list_size
1203 allocate(cs%DL(cs%list_size+1))
1205 vol = 0.0 ; area = 0.0
1206 dprev = dlist(indx2(mls))
1212 vol = vol + area * (dprev - dlist(i))
1213 area = area + arealist(i)
1215 add_to_list = .false.
1216 if ((kl == 0) .or. (k==1))
then 1217 add_to_list = .true.
1218 elseif (dlist(indx2(k-1)) < d_list_prev-cs%D_list_min_inc)
then 1219 add_to_list = .true.
1220 d_list_prev = dlist(indx2(k-1))
1223 if (add_to_list)
then 1225 cs%DL(kl)%depth = dlist(i)
1226 cs%DL(kl)%area = area
1227 cs%DL(kl)%vol_below = vol
1232 do while (kl < list_size)
1235 cs%DL(kl)%vol_below = cs%DL(kl-1)%vol_below * 1.000001
1236 cs%DL(kl)%area = cs%DL(kl-1)%area
1237 cs%DL(kl)%depth = cs%DL(kl-1)%depth
1240 cs%DL(cs%list_size+1)%vol_below = cs%DL(cs%list_size)%vol_below * 1000.0
1241 cs%DL(cs%list_size+1)%area = cs%DL(cs%list_size)%area
1242 cs%DL(cs%list_size+1)%depth = cs%DL(cs%list_size)%depth
1244 end subroutine create_depth_list
1247 subroutine write_depth_list(G, US, CS, filename, list_size)
1252 character(len=*),
intent(in) :: filename
1253 integer,
intent(in) :: list_size
1255 real,
allocatable :: tmp(:)
1256 integer :: ncid, dimid(1), Did, Aid, Vid, status, k
1257 character(len=16) :: depth_chksum, area_chksum
1260 call get_depth_list_checksums(g, depth_chksum, area_chksum)
1262 if (.not.is_root_pe())
return 1264 allocate(tmp(list_size)) ; tmp(:) = 0.0
1266 status = nf90_create(filename, 0, ncid)
1267 if (status /= nf90_noerr)
then 1268 call mom_error(warning, filename//trim(nf90_strerror(status)))
1272 status = nf90_def_dim(ncid,
"list", list_size, dimid(1))
1273 if (status /= nf90_noerr)
call mom_error(warning, &
1274 filename//trim(nf90_strerror(status)))
1276 status = nf90_def_var(ncid,
"depth", nf90_double, dimid, did)
1277 if (status /= nf90_noerr)
call mom_error(warning, &
1278 filename//
" depth "//trim(nf90_strerror(status)))
1279 status = nf90_put_att(ncid, did,
"long_name",
"Sorted depth")
1280 if (status /= nf90_noerr)
call mom_error(warning, &
1281 filename//
" depth "//trim(nf90_strerror(status)))
1282 status = nf90_put_att(ncid, did,
"units",
"m")
1283 if (status /= nf90_noerr)
call mom_error(warning, &
1284 filename//
" depth "//trim(nf90_strerror(status)))
1286 status = nf90_def_var(ncid,
"area", nf90_double, dimid, aid)
1287 if (status /= nf90_noerr)
call mom_error(warning, &
1288 filename//
" area "//trim(nf90_strerror(status)))
1289 status = nf90_put_att(ncid, aid,
"long_name",
"Open area at depth")
1290 if (status /= nf90_noerr)
call mom_error(warning, &
1291 filename//
" area "//trim(nf90_strerror(status)))
1292 status = nf90_put_att(ncid, aid,
"units",
"m2")
1293 if (status /= nf90_noerr)
call mom_error(warning, &
1294 filename//
" area "//trim(nf90_strerror(status)))
1296 status = nf90_def_var(ncid,
"vol_below", nf90_double, dimid, vid)
1297 if (status /= nf90_noerr)
call mom_error(warning, &
1298 filename//
" vol_below "//trim(nf90_strerror(status)))
1299 status = nf90_put_att(ncid, vid,
"long_name",
"Open volume below depth")
1300 if (status /= nf90_noerr)
call mom_error(warning, &
1301 filename//
" vol_below "//trim(nf90_strerror(status)))
1302 status = nf90_put_att(ncid, vid,
"units",
"m3")
1303 if (status /= nf90_noerr)
call mom_error(warning, &
1304 filename//
" vol_below "//trim(nf90_strerror(status)))
1307 status = nf90_put_att(ncid, nf90_global, depth_chksum_attr, depth_chksum)
1308 if (status /= nf90_noerr)
call mom_error(warning, &
1309 filename//
" "//depth_chksum_attr//
" "//trim(nf90_strerror(status)))
1311 status = nf90_put_att(ncid, nf90_global, area_chksum_attr, area_chksum)
1312 if (status /= nf90_noerr)
call mom_error(warning, &
1313 filename//
" "//area_chksum_attr//
" "//trim(nf90_strerror(status)))
1315 status = nf90_enddef(ncid)
1316 if (status /= nf90_noerr)
call mom_error(warning, &
1317 filename//trim(nf90_strerror(status)))
1319 do k=1,list_size ; tmp(k) = us%Z_to_m*cs%DL(k)%depth ;
enddo 1320 status = nf90_put_var(ncid, did, tmp)
1321 if (status /= nf90_noerr)
call mom_error(warning, &
1322 filename//
" depth "//trim(nf90_strerror(status)))
1324 do k=1,list_size ; tmp(k) = us%L_to_m**2*cs%DL(k)%area ;
enddo 1325 status = nf90_put_var(ncid, aid, tmp)
1326 if (status /= nf90_noerr)
call mom_error(warning, &
1327 filename//
" area "//trim(nf90_strerror(status)))
1329 do k=1,list_size ; tmp(k) = us%Z_to_m*us%L_to_m**2*cs%DL(k)%vol_below ;
enddo 1330 status = nf90_put_var(ncid, vid, tmp)
1331 if (status /= nf90_noerr)
call mom_error(warning, &
1332 filename//
" vol_below "//trim(nf90_strerror(status)))
1334 status = nf90_close(ncid)
1335 if (status /= nf90_noerr)
call mom_error(warning, &
1336 filename//trim(nf90_strerror(status)))
1338 end subroutine write_depth_list
1342 subroutine read_depth_list(G, US, CS, filename)
1347 character(len=*),
intent(in) :: filename
1349 character(len=32) :: mdl
1350 character(len=240) :: var_name, var_msg
1351 real,
allocatable :: tmp(:)
1352 integer :: ncid, status, varid, list_size, k
1353 integer :: ndim, len, var_dim_ids(nf90_max_var_dims)
1354 character(len=16) :: depth_file_chksum, depth_grid_chksum
1355 character(len=16) :: area_file_chksum, area_grid_chksum
1356 integer :: depth_attr_status, area_attr_status
1358 mdl =
"MOM_sum_output read_depth_list:" 1360 status = nf90_open(filename, nf90_nowrite, ncid)
1361 if (status /= nf90_noerr)
then 1362 call mom_error(fatal,mdl//
" Difficulties opening "//trim(filename)// &
1363 " - "//trim(nf90_strerror(status)))
1367 depth_attr_status = nf90_get_att(ncid, nf90_global, depth_chksum_attr, &
1369 area_attr_status = nf90_get_att(ncid, nf90_global, area_chksum_attr, &
1372 if (any([depth_attr_status, area_attr_status] == nf90_enotatt))
then 1373 var_msg = trim(cs%depth_list_file) //
" checksums are missing;" 1374 if (cs%require_depth_list_chksum)
then 1375 call mom_error(fatal, trim(var_msg) //
" aborting.")
1376 elseif (cs%update_depth_list_chksum)
then 1377 call mom_error(warning, trim(var_msg) //
" updating file.")
1378 call create_depth_list(g, cs)
1379 call write_depth_list(g, us, cs, cs%depth_list_file, cs%list_size+1)
1382 call mom_error(warning, &
1383 trim(var_msg) //
" some diagnostics may not be reproducible.")
1387 if (depth_attr_status /= nf90_noerr)
then 1388 var_msg = mdl //
"Failed to read " // trim(filename) //
":" &
1389 // depth_chksum_attr
1390 call mom_error(fatal, &
1391 trim(var_msg) //
" - " // nf90_strerror(depth_attr_status))
1394 if (area_attr_status /= nf90_noerr)
then 1395 var_msg = mdl //
"Failed to read " // trim(filename) //
":" &
1397 call mom_error(fatal, &
1398 trim(var_msg) //
" - " // nf90_strerror(area_attr_status))
1401 call get_depth_list_checksums(g, depth_grid_chksum, area_grid_chksum)
1403 if (depth_grid_chksum /= depth_file_chksum &
1404 .or. area_grid_chksum /= area_file_chksum)
then 1405 var_msg = trim(cs%depth_list_file) //
" checksums do not match;" 1406 if (cs%require_depth_list_chksum)
then 1407 call mom_error(fatal, trim(var_msg) //
" aborting.")
1408 elseif (cs%update_depth_list_chksum)
then 1409 call mom_error(warning, trim(var_msg) //
" updating file.")
1410 call create_depth_list(g, cs)
1411 call write_depth_list(g, us, cs, cs%depth_list_file, cs%list_size+1)
1414 call mom_error(warning, &
1415 trim(var_msg) //
" some diagnostics may not be reproducible.")
1421 var_msg = trim(var_name)//
" in "//trim(filename)//
" - " 1422 status = nf90_inq_varid(ncid, var_name, varid)
1423 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1424 " Difficulties finding variable "//trim(var_msg)//&
1425 trim(nf90_strerror(status)))
1427 status = nf90_inquire_variable(ncid, varid, ndims=ndim, dimids=var_dim_ids)
1428 if (status /= nf90_noerr)
then 1429 call mom_error(fatal,mdl//
" cannot inquire about "//trim(var_msg)//&
1430 trim(nf90_strerror(status)))
1431 elseif (ndim > 1)
then 1432 call mom_error(fatal,mdl//
" "//trim(var_msg)//&
1433 " has too many or too few dimensions.")
1437 status = nf90_inquire_dimension(ncid, var_dim_ids(1), len=list_size)
1438 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1439 " cannot inquire about dimension(1) of "//trim(var_msg)//&
1440 trim(nf90_strerror(status)))
1442 cs%list_size = list_size-1
1443 allocate(cs%DL(list_size))
1444 allocate(tmp(list_size))
1446 status = nf90_get_var(ncid, varid, tmp)
1447 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1448 " Difficulties reading variable "//trim(var_msg)//&
1449 trim(nf90_strerror(status)))
1451 do k=1,list_size ; cs%DL(k)%depth = us%m_to_Z*tmp(k) ;
enddo 1454 var_msg = trim(var_name)//
" in "//trim(filename)//
" - " 1455 status = nf90_inq_varid(ncid, var_name, varid)
1456 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1457 " Difficulties finding variable "//trim(var_msg)//&
1458 trim(nf90_strerror(status)))
1459 status = nf90_get_var(ncid, varid, tmp)
1460 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1461 " Difficulties reading variable "//trim(var_msg)//&
1462 trim(nf90_strerror(status)))
1464 do k=1,list_size ; cs%DL(k)%area = us%m_to_L**2*tmp(k) ;
enddo 1466 var_name =
"vol_below" 1467 var_msg = trim(var_name)//
" in "//trim(filename)
1468 status = nf90_inq_varid(ncid, var_name, varid)
1469 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1470 " Difficulties finding variable "//trim(var_msg)//&
1471 trim(nf90_strerror(status)))
1472 status = nf90_get_var(ncid, varid, tmp)
1473 if (status /= nf90_noerr)
call mom_error(fatal,mdl// &
1474 " Difficulties reading variable "//trim(var_msg)//&
1475 trim(nf90_strerror(status)))
1477 do k=1,list_size ; cs%DL(k)%vol_below = us%m_to_Z*us%m_to_L**2*tmp(k) ;
enddo 1479 status = nf90_close(ncid)
1480 if (status /= nf90_noerr)
call mom_error(warning, mdl// &
1481 " Difficulties closing "//trim(filename)//
" - "//trim(nf90_strerror(status)))
1485 end subroutine read_depth_list
1498 subroutine get_depth_list_checksums(G, depth_chksum, area_chksum)
1500 character(len=16),
intent(out) :: depth_chksum
1501 character(len=16),
intent(out) :: area_chksum
1504 real,
allocatable :: field(:,:)
1506 allocate(field(g%isc:g%iec, g%jsc:g%jec))
1509 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1510 field(i,j) = g%bathyT(i,j)
1512 write(depth_chksum,
'(Z16)') mpp_chksum(field(:,:))
1515 do j=g%jsc,g%jec ;
do i=g%isc,g%iec
1516 field(i,j) = g%mask2dT(i,j) * g%US%L_to_m**2*g%areaT(i,j)
1518 write(area_chksum,
'(Z16)') mpp_chksum(field(:,:))
1521 end subroutine get_depth_list_checksums
Calculates the heights of sruface or all interfaces from layer thicknesses.
Pointers to various fields which may be used describe the surface state of MOM, and which will be ret...
Wraps the FMS time manager functions.
This module implements boundary forcing for MOM6.
Reports integrated quantities for monitoring the model state.
Functions for calculating interface heights, including free surface height.
Ocean grid type. See mom_grid for details.
Sum a value or 1-d array of values across processors, returning the sums in place.
A structure that can be parsed to read and document run-time parameters.
Provides the ocean grid type.
Open boundary segment data structure.
A list of depths and corresponding globally integrated ocean area at each depth and the ocean volume ...
This module contains I/O framework code.
The MOM6 facility to parse input files for runtime parameters.
An overloaded interface to log the values of various types of parameters.
Orchestrates the registration and calling of tracer packages.
Interfaces to non-domain-oriented communication subroutines, including the MOM6 reproducing sums faci...
The control structure for the MOM_sum_output module.
Describes various unit conversion factors.
Provides a transparent unit rescaling type to facilitate dimensional consistency testing.
Routines for error handling and I/O management.
The control structure for orchestrating the calling of tracer packages.
Structure that contains pointers to the boundary forcing used to drive the liquid ocean simulated by ...
Type for describing a variable, typically a tracer.
An overloaded interface to log version information about modules.
Describes the vertical ocean grid, including unit conversion factors.
Indicate whether a file exists, perhaps with domain decomposition.
Pointers to an assortment of thermodynamic fields that may be available, including potential temperat...
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 ...
Controls where open boundary conditions are applied.
Provides a transparent vertical ocean grid type and supporting routines.
Provides transparent structures with groups of MOM6 variables and supporting routines.
An overloaded interface to read and log the values of various types of parameters.
The Extended Fixed Point (EFP) type provides a public interface for doing sums and taking differences...