forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 63
/
Copy pathMOM_surface_forcing_gfdl.F90
1850 lines (1645 loc) · 101 KB
/
MOM_surface_forcing_gfdl.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
module MOM_surface_forcing_gfdl
! This file is part of MOM6. See LICENSE.md for the license.
!#CTRL# use MOM_controlled_forcing, only : apply_ctrl_forcing, register_ctrl_forcing_restarts
!#CTRL# use MOM_controlled_forcing, only : controlled_forcing_init, controlled_forcing_end
!#CTRL# use MOM_controlled_forcing, only : ctrl_forcing_CS
use MOM_coms, only : reproducing_sum, field_chksum
use MOM_constants, only : hlv, hlf
use MOM_coupler_types, only : coupler_2d_bc_type, coupler_type_write_chksums
use MOM_coupler_types, only : coupler_type_initialized, coupler_type_spawn
use MOM_coupler_types, only : coupler_type_copy_data
use MOM_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
use MOM_cpu_clock, only : CLOCK_SUBCOMPONENT
use MOM_data_override, only : data_override_init, data_override
use MOM_diag_mediator, only : diag_ctrl, safe_alloc_ptr, time_type
use MOM_domains, only : pass_vector, pass_var, fill_symmetric_edges
use MOM_domains, only : AGRID, BGRID_NE, CGRID_NE, To_All
use MOM_domains, only : To_North, To_East, Omit_Corners
use MOM_error_handler, only : MOM_error, WARNING, FATAL, is_root_pe, MOM_mesg
use MOM_file_parser, only : get_param, log_param, log_version, param_file_type
use MOM_forcing_type, only : forcing, mech_forcing
use MOM_forcing_type, only : forcing_diags, mech_forcing_diags, register_forcing_type_diags
use MOM_forcing_type, only : allocate_forcing_type, deallocate_forcing_type
use MOM_forcing_type, only : allocate_mech_forcing, deallocate_mech_forcing
use MOM_get_input, only : Get_MOM_Input, directories
use MOM_grid, only : ocean_grid_type
use MOM_interpolate, only : init_external_field, time_interp_external
use MOM_interpolate, only : time_interp_external_init
use MOM_interpolate, only : external_field
use MOM_io, only : slasher, write_version_number, MOM_read_data
use MOM_io, only : read_netCDF_data
use MOM_io, only : stdout_if_root
use MOM_restart, only : register_restart_field, restart_init, MOM_restart_CS
use MOM_restart, only : restart_init_end, save_restart, restore_state
use MOM_string_functions, only : uppercase
use MOM_spatial_means, only : adjust_area_mean_to_zero
use MOM_unit_scaling, only : unit_scale_type
use MOM_variables, only : surface
use user_revise_forcing, only : user_alter_forcing, user_revise_forcing_init
use user_revise_forcing, only : user_revise_forcing_CS
use iso_fortran_env, only : int64
implicit none ; private
#include <MOM_memory.h>
public convert_IOB_to_fluxes, convert_IOB_to_forces
public surface_forcing_init
public ice_ocn_bnd_type_chksum
public forcing_save_restart
! A note on unit descriptions in comments: MOM6 uses units that can be rescaled for dimensional
! consistency testing. These are noted in comments with units like Z, H, L, and T, along with
! their mks counterparts with notation like "a velocity [Z T-1 ~> m s-1]". If the units
! vary with the Boussinesq approximation, the Boussinesq variant is given first.
!> surface_forcing_CS is a structure containing pointers to the forcing fields
!! which may be used to drive MOM. All fluxes are positive downward.
type, public :: surface_forcing_CS ; private
integer :: wind_stagger !< AGRID, BGRID_NE, or CGRID_NE (integer values
!! from MOM_domains) to indicate the staggering of
!! the winds that are being provided in calls to
!! update_ocean_model.
logical :: use_temperature !< If true, temp and saln used as state variables.
logical :: nonBous !< If true, this run is fully non-Boussinesq
real :: wind_stress_multiplier !< A multiplier applied to incoming wind stress [nondim].
real :: Rho0 !< Boussinesq reference density [R ~> kg m-3]
real :: area_surf = -1.0 !< Total ocean surface area [m2]
real :: latent_heat_fusion !< Latent heat of fusion [Q ~> J kg-1]
real :: latent_heat_vapor !< Latent heat of vaporization [Q ~> J kg-1]
real :: max_p_surf !< The maximum surface pressure that can be exerted by
!! the atmosphere and floating sea-ice [R L2 T-2 ~> Pa].
!! This is needed because the FMS coupling structure
!! does not limit the water that can be frozen out
!! of the ocean and the ice-ocean heat fluxes are
!! treated explicitly.
logical :: use_limited_P_SSH !< If true, return the sea surface height with
!! the correction for the atmospheric (and sea-ice)
!! pressure limited by max_p_surf instead of the
!! full atmospheric pressure. The default is true.
logical :: approx_net_mass_src !< If true, use the net mass sources from the ice-ocean boundary
!! type without any further adjustments to drive the ocean dynamics.
!! The actual net mass source may differ due to corrections.
real :: gust_const !< Constant unresolved background gustiness for ustar [R L Z T-2 ~> Pa]
logical :: read_gust_2d !< If true, use a 2-dimensional gustiness supplied from an input file.
real, pointer, dimension(:,:) :: &
TKE_tidal => NULL() !< Turbulent kinetic energy introduced to the bottom boundary layer
!! by drag on the tidal flows [R Z3 T-3 ~> W m-2].
real, pointer, dimension(:,:) :: &
gust => NULL() !< A spatially varying unresolved background gustiness that
!! contributes to ustar [R L Z T-2 ~> Pa]. gust is used when read_gust_2d is true.
real, pointer, dimension(:,:) :: &
ustar_tidal => NULL() !< Tidal contribution to the bottom friction velocity [Z T-1 ~> m s-1]
real :: cd_tides !< Drag coefficient that applies to the tides [nondim]
real :: utide !< Constant tidal velocity to use if read_tideamp is false [Z T-1 ~> m s-1].
logical :: read_tideamp !< If true, spatially varying tidal amplitude read from a file.
logical :: rigid_sea_ice !< If true, sea-ice exerts a rigidity that acts to damp surface
!! deflections (especially surface gravity waves). The default is false.
real :: g_Earth !< Gravitational acceleration [L2 Z-1 T-2 ~> m s-2]
real :: Kv_sea_ice !< Viscosity in sea-ice that resists sheared vertical motions [L4 Z-2 T-1 ~> m2 s-1]
real :: density_sea_ice !< Typical density of sea-ice [R ~> kg m-3]. The value is only used to convert
!! the ice pressure into appropriate units for use with Kv_sea_ice.
real :: rigid_sea_ice_mass !< A mass per unit area of sea-ice beyond which sea-ice viscosity
!! becomes effective [R Z ~> kg m-2], typically of order 1000 kg m-2.
logical :: allow_flux_adjustments !< If true, use data_override to obtain flux adjustments
logical :: restore_salt !< If true, the coupled MOM driver adds a term to restore surface
!! salinity to a specified value.
logical :: restore_temp !< If true, the coupled MOM driver adds a term to restore sea
!! surface temperature to a specified value.
real :: Flux_const_salt !< Piston velocity for surface salinity restoring [Z T-1 ~> m s-1]
real :: Flux_const_temp !< Piston velocity for surface temperature restoring [Z T-1 ~> m s-1]
real :: rho_restore !< The density that is used to convert piston velocities into salt
!! or heat fluxes with salinity or temperature restoring [R ~> kg m-3]
logical :: trestore_SPEAR_ECDA !< If true, modify restoring data wrt local SSS
real :: SPEAR_dTf_dS !< The derivative of the freezing temperature with
!! salinity [C S-1 ~> degC ppt-1].
logical :: salt_restore_as_sflux !< If true, SSS restore as salt flux instead of water flux
logical :: adjust_net_srestore_to_zero !< Adjust srestore to zero (for both salt_flux or vprec)
logical :: adjust_net_srestore_by_scaling !< Adjust srestore w/o moving zero contour
logical :: adjust_net_fresh_water_to_zero !< Adjust net surface fresh-water (with restoring) to zero
logical :: use_net_FW_adjustment_sign_bug !< Use the wrong sign when adjusting net FW
logical :: adjust_net_fresh_water_by_scaling !< Adjust net surface fresh-water w/o moving zero contour
logical :: mask_srestore_under_ice !< If true, use an ice mask defined by frazil criteria
!! for salinity restoring.
real :: ice_salt_concentration !< Salt concentration for sea ice [kg/kg]
logical :: mask_srestore_marginal_seas !< If true, then mask SSS restoring in marginal seas
real :: max_delta_srestore !< Maximum delta salinity used for restoring [S ~> ppt]
real :: max_delta_trestore !< Maximum delta sst used for restoring [C ~> degC]
real, pointer, dimension(:,:) :: basin_mask => NULL() !< Mask for surface salinity restoring by basin [nondim]
integer :: answer_date !< The vintage of the order of arithmetic and expressions in the
!! gustiness calculations. Values below 20190101 recover the answers
!! from the end of 2018, while higher values use a simpler expression
!! to calculate gustiness.
logical :: ustar_gustless_bug !< If true, include a bug in the time-averaging of the
!! gustless wind friction velocity.
logical :: check_no_land_fluxes !< Return warning if IOB flux over land is non-zero
type(diag_ctrl), pointer :: diag => NULL() !< Structure to regulate diagnostic output timing
character(len=200) :: inputdir !< Directory where NetCDF input files are
character(len=200) :: salt_restore_file !< Filename for salt restoring data
character(len=30) :: salt_restore_var_name !< Name of surface salinity in salt_restore_file
logical :: mask_srestore !< If true, apply a 2-dimensional mask to the surface
!! salinity restoring fluxes. The masking file should be
!! in inputdir/salt_restore_mask.nc and the field should
!! be named 'mask'
real, pointer, dimension(:,:) :: srestore_mask => NULL() !< mask for SSS restoring [nondim]
character(len=200) :: temp_restore_file !< Filename for sst restoring data
character(len=30) :: temp_restore_var_name !< Name of surface temperature in temp_restore_file
logical :: mask_trestore !< If true, apply a 2-dimensional mask to the surface
!! temperature restoring fluxes. The masking file should be
!! in inputdir/temp_restore_mask.nc and the field should
!! be named 'mask'
real, pointer, dimension(:,:) :: trestore_mask => NULL() !< Mask for SST restoring [nondim]
type(external_field) :: srestore_handle
!< Handle for time-interpolated salt restoration field
type(external_field) :: trestore_handle
!< Handle for time-interpolated temperature restoration field
type(forcing_diags), public :: handles !< Diagnostics handles
!#CTRL# type(ctrl_forcing_CS), pointer :: ctrl_forcing_CSp => NULL()
type(MOM_restart_CS), pointer :: restart_CSp => NULL() !< A pointer to the restart control structure
type(user_revise_forcing_CS), pointer :: urf_CS => NULL() !< A control structure for user forcing revisions
end type surface_forcing_CS
!> ice_ocean_boundary_type is a structure corresponding to forcing, but with the elements, units,
!! and conventions that exactly conform to the use for MOM6-based coupled models.
type, public :: ice_ocean_boundary_type
real, pointer, dimension(:,:) :: u_flux =>NULL() !< i-direction wind stress [Pa]
real, pointer, dimension(:,:) :: v_flux =>NULL() !< j-direction wind stress [Pa]
real, pointer, dimension(:,:) :: t_flux =>NULL() !< sensible heat flux [W m-2]
real, pointer, dimension(:,:) :: q_flux =>NULL() !< specific humidity flux [kg m-2 s-1]
real, pointer, dimension(:,:) :: salt_flux =>NULL() !< salt flux [kg m-2 s-1]
real, pointer, dimension(:,:) :: excess_salt =>NULL() !< salt left behind by brine rejection [kg m-2 s-1]
real, pointer, dimension(:,:) :: lw_flux =>NULL() !< long wave radiation [W m-2]
real, pointer, dimension(:,:) :: sw_flux_vis_dir =>NULL() !< direct visible sw radiation [W m-2]
real, pointer, dimension(:,:) :: sw_flux_vis_dif =>NULL() !< diffuse visible sw radiation [W m-2]
real, pointer, dimension(:,:) :: sw_flux_nir_dir =>NULL() !< direct Near InfraRed sw radiation [W m-2]
real, pointer, dimension(:,:) :: sw_flux_nir_dif =>NULL() !< diffuse Near InfraRed sw radiation [W m-2]
real, pointer, dimension(:,:) :: lprec =>NULL() !< mass flux of liquid precip [kg m-2 s-1]
real, pointer, dimension(:,:) :: fprec =>NULL() !< mass flux of frozen precip [kg m-2 s-1]
real, pointer, dimension(:,:) :: seaice_melt =>NULL() !< mass flux of melted sea ice [kg m-2 s-1]
real, pointer, dimension(:,:) :: runoff =>NULL() !< mass flux of liquid runoff [kg m-2 s-1]
real, pointer, dimension(:,:) :: calving =>NULL() !< mass flux of frozen runoff [kg m-2 s-1]
real, pointer, dimension(:,:) :: stress_mag =>NULL() !< The time-mean magnitude of the stress on the ocean [Pa]
real, pointer, dimension(:,:) :: ustar_berg =>NULL() !< frictional velocity beneath icebergs [m s-1]
real, pointer, dimension(:,:) :: area_berg =>NULL() !< fractional area covered by icebergs [m2 m-2]
real, pointer, dimension(:,:) :: mass_berg =>NULL() !< mass of icebergs per unit ocean area [kg m-2]
real, pointer, dimension(:,:) :: runoff_hflx =>NULL() !< heat content of liquid runoff [W m-2]
real, pointer, dimension(:,:) :: calving_hflx =>NULL() !< heat content of frozen runoff [W m-2]
real, pointer, dimension(:,:) :: p =>NULL() !< pressure of overlying ice and atmosphere
!< on ocean surface [Pa]
real, pointer, dimension(:,:) :: mi =>NULL() !< mass of ice per unit ocean area [kg m-2]
real, pointer, dimension(:,:) :: ice_rigidity =>NULL() !< rigidity of the sea ice, sea-ice and
!! ice-shelves, expressed as a coefficient
!! for divergence damping, as determined
!! outside of the ocean model [m3 s-1]
integer :: xtype !< The type of the exchange - REGRID, REDIST or DIRECT
type(coupler_2d_bc_type) :: fluxes !< A structure that may contain an array of named fields
!! used for passive tracer fluxes.
integer :: wind_stagger = -999 !< A flag indicating the spatial discretization of wind stresses.
!! This flag may be set by the flux-exchange code, based on what
!! the sea-ice model is providing. Otherwise, the value from
!! the surface_forcing_CS is used.
end type ice_ocean_boundary_type
integer :: id_clock_forcing !< A CPU time clock
contains
!> This subroutine translates the Ice_ocean_boundary_type into a MOM
!! thermodynamic forcing type, including changes of units, sign conventions,
!! and putting the fields into arrays with MOM-standard halos.
subroutine convert_IOB_to_fluxes(IOB, fluxes, index_bounds, Time, valid_time, G, US, CS, sfc_state)
type(ice_ocean_boundary_type), &
target, intent(in) :: IOB !< An ice-ocean boundary type with fluxes to drive
!! the ocean in a coupled model
type(forcing), intent(inout) :: fluxes !< A structure containing pointers to all
!! possible mass, heat or salt flux forcing fields.
!! Unused fields have NULL ptrs.
integer, dimension(4), intent(in) :: index_bounds !< The i- and j- size of the arrays in IOB.
type(time_type), intent(in) :: Time !< The time of the fluxes, used for interpolating the
!! salinity to the right time, when it is being restored.
real, intent(in) :: valid_time !< The amount of time over which these fluxes
!! should be applied [T ~> s].
type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(surface_forcing_CS),pointer :: CS !< A pointer to the control structure returned by a
!! previous call to surface_forcing_init.
type(surface), intent(in) :: sfc_state !< A structure containing fields that describe the
!! surface state of the ocean.
real, dimension(SZI_(G),SZJ_(G)) :: &
data_restore, & ! The surface value toward which to restore [S ~> ppt] or [C ~> degC]
SST_anom, & ! Instantaneous sea surface temperature anomalies from a target value [C ~> degC]
SSS_anom, & ! Instantaneous sea surface salinity anomalies from a target value [S ~> ppt]
SSS_mean, & ! A (mean?) salinity about which to normalize local salinity
! anomalies when calculating restorative precipitation anomalies [S ~> ppt]
net_FW, & ! The area integrated net freshwater flux into the ocean [kg s-1]
net_FW2, & ! The net freshwater flux into the ocean [kg m-2 s-1]
work_sum, & ! A 2-d array that is used as the work space for global sums [m2] or [kg s-1]
open_ocn_mask ! a binary field indicating where ice is present based on frazil criteria [nondim]
integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq, i0, j0
integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, isr, ier, jsr, jer
integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
real :: delta_sss ! temporary storage for sss diff from restoring value [S ~> ppt]
real :: delta_sst ! temporary storage for sst diff from restoring value [C ~> degC]
real :: kg_m2_s_conversion ! A combination of unit conversion factors for rescaling
! mass fluxes [R Z s m2 kg-1 T-1 ~> 1]
real :: rhoXcp ! Reference density times heat capacity times unit scaling
! factors [Q R C-1 ~> J m-3 degC-1]
real :: sign_for_net_FW_bug ! Should be +1. but an old bug can be recovered by using -1 [nondim]
call cpu_clock_begin(id_clock_forcing)
isc_bnd = index_bounds(1) ; iec_bnd = index_bounds(2)
jsc_bnd = index_bounds(3) ; jec_bnd = index_bounds(4)
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
Isq = G%IscB ; Ieq = G%IecB ; Jsq = G%JscB ; Jeq = G%JecB
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
IsdB = G%IsdB ; IedB = G%IedB ; JsdB = G%JsdB ; JedB = G%JedB
isr = is-isd+1 ; ier = ie-isd+1 ; jsr = js-jsd+1 ; jer = je-jsd+1
kg_m2_s_conversion = US%kg_m2s_to_RZ_T
if (CS%restore_temp) rhoXcp = CS%rho_restore * fluxes%C_p
open_ocn_mask(:,:) = 1.0
fluxes%vPrecGlobalAdj = 0.0
fluxes%vPrecGlobalScl = 0.0
fluxes%saltFluxGlobalAdj = 0.0
fluxes%saltFluxGlobalScl = 0.0
fluxes%netFWGlobalAdj = 0.0
fluxes%netFWGlobalScl = 0.0
! allocation and initialization if this is the first time that this
! flux type has been used.
if (fluxes%dt_buoy_accum < 0) then
call allocate_forcing_type(G, fluxes, water=.true., heat=.true., ustar=.not.CS%nonBous, press=.true., &
fix_accum_bug=.not.CS%ustar_gustless_bug, tau_mag=CS%nonBous)
call safe_alloc_ptr(fluxes%sw_vis_dir,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%sw_vis_dif,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%sw_nir_dir,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%sw_nir_dif,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%p_surf,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%p_surf_full,isd,ied,jsd,jed)
if (CS%use_limited_P_SSH) then
fluxes%p_surf_SSH => fluxes%p_surf
else
fluxes%p_surf_SSH => fluxes%p_surf_full
endif
call safe_alloc_ptr(fluxes%salt_flux,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%salt_flux_in,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%TKE_tidal,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%ustar_tidal,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%heat_added,isd,ied,jsd,jed)
call safe_alloc_ptr(fluxes%salt_flux_added,isd,ied,jsd,jed)
if (associated(IOB%excess_salt)) call safe_alloc_ptr(fluxes%salt_left_behind,isd,ied,jsd,jed)
do j=js-2,je+2 ; do i=is-2,ie+2
fluxes%TKE_tidal(i,j) = CS%TKE_tidal(i,j)
fluxes%ustar_tidal(i,j) = CS%ustar_tidal(i,j)
enddo ; enddo
endif ! endif for allocation and initialization
if (((associated(IOB%ustar_berg) .and. (.not.associated(fluxes%ustar_berg))) &
.or. (associated(IOB%area_berg) .and. (.not.associated(fluxes%area_berg)))) &
.or. (associated(IOB%mass_berg) .and. (.not.associated(fluxes%mass_berg)))) &
call allocate_forcing_type(G, fluxes, iceberg=.true.)
if ((.not.coupler_type_initialized(fluxes%tr_fluxes)) .and. &
coupler_type_initialized(IOB%fluxes)) &
call coupler_type_spawn(IOB%fluxes, fluxes%tr_fluxes, (/is,is,ie,ie/), (/js,js,je,je/))
! It might prove valuable to use the same array extents as the rest of the
! ocean model, rather than using haloless arrays, in which case the last line
! would be: ( (/isd,is,ie,ied/), (/jsd,js,je,jed/))
! allocation and initialization on first call to this routine
if (CS%area_surf < 0.0) then
do j=js,je ; do i=is,ie
work_sum(i,j) = US%L_to_m**2*G%areaT(i,j) * G%mask2dT(i,j)
enddo ; enddo
CS%area_surf = reproducing_sum(work_sum, isr, ier, jsr, jer)
endif ! endif for allocation and initialization
! Indicate that there are new unused fluxes.
fluxes%fluxes_used = .false.
fluxes%dt_buoy_accum = valid_time
fluxes%heat_added(:,:) = 0.0
fluxes%salt_flux_added(:,:) = 0.0
do j=js,je ; do i=is,ie
fluxes%salt_flux(i,j) = 0.0
fluxes%vprec(i,j) = 0.0
enddo ; enddo
! Salinity restoring logic
if (CS%restore_salt) then
call time_interp_external(CS%srestore_handle, Time, data_restore, scale=US%ppt_to_S)
! open_ocn_mask indicates where to restore salinity (1 means restore, 0 does not)
open_ocn_mask(:,:) = 1.0
if (CS%mask_srestore_under_ice) then ! Do not restore under sea-ice
do j=js,je ; do i=is,ie
if (sfc_state%SST(i,j) <= CS%SPEAR_dTf_dS*sfc_state%SSS(i,j)) open_ocn_mask(i,j)=0.0
enddo ; enddo
endif
if (CS%salt_restore_as_sflux) then
do j=js,je ; do i=is,ie
delta_sss = data_restore(i,j) - sfc_state%SSS(i,j)
delta_sss = sign(1.0,delta_sss) * min(abs(delta_sss), CS%max_delta_srestore)
fluxes%salt_flux(i,j) = 1.e-3*US%S_to_ppt*G%mask2dT(i,j) * (CS%rho_restore*CS%Flux_const_salt)* &
(CS%basin_mask(i,j)*open_ocn_mask(i,j)*CS%srestore_mask(i,j)) * delta_sss ! R Z T-1 ~> kg Salt m-2 s-1
enddo ; enddo
if (CS%adjust_net_srestore_to_zero) then
if (CS%adjust_net_srestore_by_scaling) then
call adjust_area_mean_to_zero(fluxes%salt_flux, G, fluxes%saltFluxGlobalScl, &
unit_scale=US%RZ_T_to_kg_m2s)
fluxes%saltFluxGlobalAdj = 0.
else
work_sum(is:ie,js:je) = US%L_to_m**2*US%RZ_T_to_kg_m2s * &
G%areaT(is:ie,js:je)*fluxes%salt_flux(is:ie,js:je) * G%mask2dT(is:ie,js:je)
fluxes%saltFluxGlobalAdj = reproducing_sum(work_sum(:,:), isr,ier, jsr,jer)/CS%area_surf
fluxes%salt_flux(is:ie,js:je) = fluxes%salt_flux(is:ie,js:je) - &
kg_m2_s_conversion * fluxes%saltFluxGlobalAdj * G%mask2dT(is:ie,js:je)
endif
endif
fluxes%salt_flux_added(is:ie,js:je) = fluxes%salt_flux(is:ie,js:je) ! Diagnostic
else
do j=js,je ; do i=is,ie
if (G%mask2dT(i,j) > 0.0) then
delta_sss = sfc_state%SSS(i,j) - data_restore(i,j)
delta_sss = sign(1.0,delta_sss) * min(abs(delta_sss), CS%max_delta_srestore)
fluxes%vprec(i,j) = (CS%basin_mask(i,j)*open_ocn_mask(i,j)*CS%srestore_mask(i,j))* &
(CS%rho_restore*CS%Flux_const_salt) * &
delta_sss / (0.5*(sfc_state%SSS(i,j) + data_restore(i,j)))
endif
enddo ; enddo
if (CS%adjust_net_srestore_to_zero) then
if (CS%adjust_net_srestore_by_scaling) then
call adjust_area_mean_to_zero(fluxes%vprec, G, fluxes%vPrecGlobalScl, &
unit_scale=US%RZ_T_to_kg_m2s)
fluxes%vPrecGlobalAdj = 0.
else
work_sum(is:ie,js:je) = US%L_to_m**2*G%areaT(is:ie,js:je) * &
US%RZ_T_to_kg_m2s*fluxes%vprec(is:ie,js:je)
fluxes%vPrecGlobalAdj = reproducing_sum(work_sum(:,:), isr, ier, jsr, jer) / CS%area_surf
do j=js,je ; do i=is,ie
fluxes%vprec(i,j) = ( fluxes%vprec(i,j) - kg_m2_s_conversion*fluxes%vPrecGlobalAdj ) * G%mask2dT(i,j)
enddo ; enddo
endif
endif
endif
endif
! SST restoring logic
if (CS%restore_temp) then
call time_interp_external(CS%trestore_handle, Time, data_restore, scale=US%degC_to_C)
if ( CS%trestore_SPEAR_ECDA ) then
do j=js,je ; do i=is,ie
if (abs(data_restore(i,j)+1.8*US%degC_to_C) < 0.0001*US%degC_to_C) then
data_restore(i,j) = CS%SPEAR_dTf_dS*sfc_state%SSS(i,j)
endif
enddo ; enddo
endif
do j=js,je ; do i=is,ie
delta_sst = data_restore(i,j) - sfc_state%SST(i,j)
delta_sst = sign(1.0,delta_sst) * min(abs(delta_sst), CS%max_delta_trestore)
fluxes%heat_added(i,j) = G%mask2dT(i,j) * CS%trestore_mask(i,j) * &
rhoXcp * delta_sst * CS%Flux_const_temp ! [Q R Z T-1 ~> W m-2]
enddo ; enddo
endif
! obtain fluxes from IOB; note the staggering of indices
i0 = is - isc_bnd ; j0 = js - jsc_bnd
do j=js,je ; do i=is,ie
if (associated(IOB%lprec)) then
fluxes%lprec(i,j) = kg_m2_s_conversion * IOB%lprec(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%lprec(i-i0,j-j0), G%mask2dT(i,j), i, j, 'lprec', G)
endif
if (associated(IOB%fprec)) then
fluxes%fprec(i,j) = kg_m2_s_conversion * IOB%fprec(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%fprec(i-i0,j-j0), G%mask2dT(i,j), i, j, 'fprec', G)
endif
if (associated(IOB%seaice_melt)) then
fluxes%seaice_melt(i,j) = kg_m2_s_conversion * IOB%seaice_melt(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%seaice_melt(i-i0,j-j0), G%mask2dT(i,j), i, j, 'seaice_melt', G)
endif
if (associated(IOB%q_flux)) then
fluxes%evap(i,j) = - kg_m2_s_conversion * IOB%q_flux(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%q_flux(i-i0,j-j0), G%mask2dT(i,j), i, j, 'q_flux', G)
endif
if (associated(IOB%runoff)) then
fluxes%lrunoff(i,j) = kg_m2_s_conversion * IOB%runoff(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%runoff(i-i0,j-j0), G%mask2dT(i,j), i, j, 'runoff', G)
endif
if (associated(IOB%calving)) then
fluxes%frunoff(i,j) = kg_m2_s_conversion * IOB%calving(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%calving(i-i0,j-j0), G%mask2dT(i,j), i, j, 'calving', G)
endif
if (associated(IOB%ustar_berg)) then
fluxes%ustar_berg(i,j) = US%m_to_Z*US%T_to_s * IOB%ustar_berg(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%ustar_berg(i-i0,j-j0), G%mask2dT(i,j), i, j, 'ustar_berg', G)
endif
if (associated(IOB%area_berg)) then
fluxes%area_berg(i,j) = IOB%area_berg(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%area_berg(i-i0,j-j0), G%mask2dT(i,j), i, j, 'area_berg', G)
endif
if (associated(IOB%mass_berg)) then
fluxes%mass_berg(i,j) = US%m_to_Z*US%kg_m3_to_R * IOB%mass_berg(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%mass_berg(i-i0,j-j0), G%mask2dT(i,j), i, j, 'mass_berg', G)
endif
if (associated(IOB%runoff_hflx)) then
fluxes%heat_content_lrunoff(i,j) = US%W_m2_to_QRZ_T * IOB%runoff_hflx(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%runoff_hflx(i-i0,j-j0), G%mask2dT(i,j), i, j, 'runoff_hflx', G)
endif
if (associated(IOB%calving_hflx)) then
fluxes%heat_content_frunoff(i,j) = US%W_m2_to_QRZ_T * IOB%calving_hflx(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%calving_hflx(i-i0,j-j0), G%mask2dT(i,j), i, j, 'calving_hflx', G)
endif
if (associated(IOB%lw_flux)) then
fluxes%LW(i,j) = US%W_m2_to_QRZ_T * IOB%lw_flux(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%lw_flux(i-i0,j-j0), G%mask2dT(i,j), i, j, 'lw_flux', G)
endif
if (associated(IOB%t_flux)) then
fluxes%sens(i,j) = -US%W_m2_to_QRZ_T* IOB%t_flux(i-i0,j-j0) * G%mask2dT(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%t_flux(i-i0,j-j0), G%mask2dT(i,j), i, j, 't_flux', G)
endif
fluxes%latent(i,j) = 0.0
if (associated(IOB%fprec)) then
fluxes%latent(i,j) = fluxes%latent(i,j) - IOB%fprec(i-i0,j-j0)*kg_m2_s_conversion * CS%latent_heat_fusion
fluxes%latent_fprec_diag(i,j) = -G%mask2dT(i,j) * IOB%fprec(i-i0,j-j0)*kg_m2_s_conversion * CS%latent_heat_fusion
endif
if (associated(IOB%calving)) then
fluxes%latent(i,j) = fluxes%latent(i,j) - IOB%calving(i-i0,j-j0)*kg_m2_s_conversion * CS%latent_heat_fusion
fluxes%latent_frunoff_diag(i,j) = -G%mask2dT(i,j) * IOB%calving(i-i0,j-j0)*kg_m2_s_conversion * &
CS%latent_heat_fusion
endif
if (associated(IOB%q_flux)) then
fluxes%latent(i,j) = fluxes%latent(i,j) - IOB%q_flux(i-i0,j-j0)*kg_m2_s_conversion * CS%latent_heat_vapor
fluxes%latent_evap_diag(i,j) = -G%mask2dT(i,j) * IOB%q_flux(i-i0,j-j0)*kg_m2_s_conversion * CS%latent_heat_vapor
endif
fluxes%latent(i,j) = G%mask2dT(i,j) * fluxes%latent(i,j)
if (associated(IOB%sw_flux_vis_dir)) then
fluxes%sw_vis_dir(i,j) = G%mask2dT(i,j) * US%W_m2_to_QRZ_T * IOB%sw_flux_vis_dir(i-i0,j-j0)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%sw_flux_vis_dir(i-i0,j-j0), G%mask2dT(i,j), i, j, 'sw_flux_vis_dir', G)
endif
if (associated(IOB%sw_flux_vis_dif)) then
fluxes%sw_vis_dif(i,j) = G%mask2dT(i,j) * US%W_m2_to_QRZ_T * IOB%sw_flux_vis_dif(i-i0,j-j0)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%sw_flux_vis_dif(i-i0,j-j0), G%mask2dT(i,j), i, j, 'sw_flux_vis_dif', G)
endif
if (associated(IOB%sw_flux_nir_dir)) then
fluxes%sw_nir_dir(i,j) = G%mask2dT(i,j) * US%W_m2_to_QRZ_T * IOB%sw_flux_nir_dir(i-i0,j-j0)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%sw_flux_nir_dir(i-i0,j-j0), G%mask2dT(i,j), i, j, 'sw_flux_nir_dir', G)
endif
if (associated(IOB%sw_flux_nir_dif)) then
fluxes%sw_nir_dif(i,j) = G%mask2dT(i,j) * US%W_m2_to_QRZ_T * IOB%sw_flux_nir_dif(i-i0,j-j0)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%sw_flux_nir_dif(i-i0,j-j0), G%mask2dT(i,j), i, j, 'sw_flux_nir_dif', G)
endif
if (CS%answer_date < 20190101) then
fluxes%sw(i,j) = fluxes%sw_vis_dir(i,j) + fluxes%sw_vis_dif(i,j) + &
fluxes%sw_nir_dir(i,j) + fluxes%sw_nir_dif(i,j)
else
fluxes%sw(i,j) = (fluxes%sw_vis_dir(i,j) + fluxes%sw_vis_dif(i,j)) + &
(fluxes%sw_nir_dir(i,j) + fluxes%sw_nir_dif(i,j))
endif
enddo ; enddo
! applied surface pressure from atmosphere and cryosphere
if (associated(IOB%p)) then
if (CS%max_p_surf >= 0.0) then
do j=js,je ; do i=is,ie
fluxes%p_surf_full(i,j) = G%mask2dT(i,j) * US%Pa_to_RL2_T2*IOB%p(i-i0,j-j0)
fluxes%p_surf(i,j) = MIN(fluxes%p_surf_full(i,j),CS%max_p_surf)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%p(i-i0,j-j0), G%mask2dT(i,j), i, j, 'p', G)
enddo ; enddo
else
do j=js,je ; do i=is,ie
fluxes%p_surf_full(i,j) = G%mask2dT(i,j) * US%Pa_to_RL2_T2*IOB%p(i-i0,j-j0)
fluxes%p_surf(i,j) = fluxes%p_surf_full(i,j)
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%p(i-i0,j-j0), G%mask2dT(i,j), i, j, 'p', G)
enddo ; enddo
endif
fluxes%accumulate_p_surf = .true. ! Multiple components may contribute to surface pressure.
endif
! more salt restoring logic
if (associated(IOB%salt_flux)) then
do j=js,je ; do i=is,ie
fluxes%salt_flux(i,j) = G%mask2dT(i,j)*(fluxes%salt_flux(i,j) - kg_m2_s_conversion*IOB%salt_flux(i-i0,j-j0))
fluxes%salt_flux_in(i,j) = G%mask2dT(i,j)*( -kg_m2_s_conversion*IOB%salt_flux(i-i0,j-j0) )
if (CS%check_no_land_fluxes) &
call check_mask_val_consistency(IOB%salt_flux(i-i0,j-j0), G%mask2dT(i,j), i, j, 'salt_flux', G)
enddo ; enddo
endif
if (associated(IOB%excess_salt)) then
do j=js,je ; do i=is,ie
fluxes%salt_left_behind(i,j) = G%mask2dT(i,j)*(kg_m2_s_conversion*IOB%excess_salt(i-i0,j-j0))
enddo ; enddo
endif
!#CTRL# if (associated(CS%ctrl_forcing_CSp)) then
!#CTRL# do j=js,je ; do i=is,ie
!#CTRL# SST_anom(i,j) = sfc_state%SST(i,j) - CS%T_Restore(i,j)
!#CTRL# SSS_anom(i,j) = sfc_state%SSS(i,j) - CS%S_Restore(i,j)
!#CTRL# SSS_mean(i,j) = 0.5*(sfc_state%SSS(i,j) + CS%S_Restore(i,j))
!#CTRL# enddo ; enddo
!#CTRL# call apply_ctrl_forcing(SST_anom, SSS_anom, SSS_mean, fluxes%heat_added, &
!#CTRL# fluxes%vprec, day, valid_time, G, US, CS%ctrl_forcing_CSp)
!#CTRL# endif
! adjust the NET fresh-water flux to zero, if flagged
if (CS%adjust_net_fresh_water_to_zero) then
sign_for_net_FW_bug = 1.
if (CS%use_net_FW_adjustment_sign_bug) sign_for_net_FW_bug = -1.
do j=js,je ; do i=is,ie
net_FW(i,j) = US%RZ_T_to_kg_m2s* &
(((fluxes%lprec(i,j) + fluxes%fprec(i,j) + fluxes%seaice_melt(i,j)) + &
(fluxes%lrunoff(i,j) + fluxes%frunoff(i,j))) + &
(fluxes%evap(i,j) + fluxes%vprec(i,j)) ) * US%L_to_m**2*G%areaT(i,j)
! The following contribution appears to be calculating the volume flux of sea-ice
! melt. This calculation is clearly WRONG if either sea-ice has variable
! salinity or the sea-ice is completely fresh.
! Bob thinks this is trying ensure the net fresh-water of the ocean + sea-ice system
! is constant.
! To do this correctly we will need a sea-ice melt field added to IOB. -AJA
if (associated(IOB%salt_flux) .and. (CS%ice_salt_concentration>0.0)) &
net_FW(i,j) = net_FW(i,j) + sign_for_net_FW_bug * US%L_to_m**2*G%areaT(i,j) * &
(IOB%salt_flux(i-i0,j-j0) / CS%ice_salt_concentration)
net_FW2(i,j) = net_FW(i,j) / (US%L_to_m**2*G%areaT(i,j))
enddo ; enddo
if (CS%adjust_net_fresh_water_by_scaling) then
call adjust_area_mean_to_zero(net_FW2, G, fluxes%netFWGlobalScl)
do j=js,je ; do i=is,ie
fluxes%vprec(i,j) = fluxes%vprec(i,j) + kg_m2_s_conversion * &
(net_FW2(i,j) - net_FW(i,j)/(US%L_to_m**2*G%areaT(i,j))) * G%mask2dT(i,j)
enddo ; enddo
else
fluxes%netFWGlobalAdj = reproducing_sum(net_FW(:,:), isr, ier, jsr, jer) / CS%area_surf
do j=js,je ; do i=is,ie
fluxes%vprec(i,j) = ( fluxes%vprec(i,j) - kg_m2_s_conversion * fluxes%netFWGlobalAdj ) * G%mask2dT(i,j)
enddo ; enddo
endif
endif
! Set the wind stresses and ustar.
if (associated(fluxes%ustar) .and. associated(fluxes%ustar_gustless) .and. associated(fluxes%tau_mag) &
.and. associated(fluxes%tau_mag_gustless) ) then
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, ustar=fluxes%ustar, &
mag_tau=fluxes%tau_mag, gustless_ustar=fluxes%ustar_gustless, &
gustless_mag_tau=fluxes%tau_mag_gustless)
else
if (associated(fluxes%ustar)) &
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, ustar=fluxes%ustar)
if (associated(fluxes%ustar_gustless)) &
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, gustless_ustar=fluxes%ustar_gustless)
if (associated(fluxes%tau_mag)) &
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, mag_tau=fluxes%tau_mag)
if (associated(fluxes%tau_mag_gustless)) &
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, gustless_mag_tau=fluxes%tau_mag_gustless)
endif
if (coupler_type_initialized(fluxes%tr_fluxes) .and. &
coupler_type_initialized(IOB%fluxes)) &
call coupler_type_copy_data(IOB%fluxes, fluxes%tr_fluxes)
if (CS%allow_flux_adjustments) then
! Apply adjustments to fluxes
call apply_flux_adjustments(G, US, CS, Time, fluxes)
endif
! Allow for user-written code to alter fluxes after all the above
call user_alter_forcing(sfc_state, fluxes, Time, G, CS%urf_CS)
call cpu_clock_end(id_clock_forcing)
end subroutine convert_IOB_to_fluxes
!> This subroutine translates the Ice_ocean_boundary_type into a MOM
!! mechanical forcing type, including changes of units, sign conventions,
!! and putting the fields into arrays with MOM-standard halos.
subroutine convert_IOB_to_forces(IOB, forces, index_bounds, Time, G, US, CS, dt_forcing, reset_avg)
type(ice_ocean_boundary_type), &
target, intent(in) :: IOB !< An ice-ocean boundary type with fluxes to drive
!! the ocean in a coupled model
type(mech_forcing), intent(inout) :: forces !< A structure with the driving mechanical forces
integer, dimension(4), intent(in) :: index_bounds !< The i- and j- size of the arrays in IOB.
type(time_type), intent(in) :: Time !< The time of the fluxes, used for interpolating the
!! salinity to the right time, when it is being restored.
type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(surface_forcing_CS),pointer :: CS !< A pointer to the control structure returned by a
!! previous call to surface_forcing_init.
real, optional, intent(in) :: dt_forcing !< A time interval over which to apply the
!! current value of ustar as a weighted running
!! average [T ~> s], or if 0 do not average ustar.
!! Missing is equivalent to 0.
logical, optional, intent(in) :: reset_avg !< If true, reset the time average.
! Local variables
real, dimension(SZI_(G),SZJ_(G)) :: &
rigidity_at_h, & ! Ice rigidity at tracer points [L4 Z-1 T-1 ~> m3 s-1]
net_mass_src, & ! A temporary of net mass sources [R Z T-1 ~> kg m-2 s-1].
ustar_tmp, & ! A temporary array of ustar values [Z T-1 ~> m s-1].
tau_mag_tmp ! A temporary array of surface stress magnitudes [R Z L T-2 ~> Pa]
real :: I_GEarth ! The inverse of the gravitational acceleration [T2 Z L-2 ~> s2 m-1]
real :: Kv_rho_ice ! (CS%Kv_sea_ice / CS%density_sea_ice) [L4 Z-2 T-1 R-1 ~> m5 s-1 kg-1]
real :: mass_ice ! mass of sea ice at a face [R Z ~> kg m-2]
real :: mass_eff ! effective mass of sea ice for rigidity [R Z ~> kg m-2]
real :: wt1, wt2 ! Relative weights of previous and current values of ustar [nondim].
real :: kg_m2_s_conversion ! A combination of unit conversion factors for rescaling
! mass fluxes [R Z s m2 kg-1 T-1 ~> 1]
integer :: i, j, is, ie, js, je, Isq, Ieq, Jsq, Jeq, i0, j0
integer :: isd, ied, jsd, jed, IsdB, IedB, JsdB, JedB, isr, ier, jsr, jer
integer :: isc_bnd, iec_bnd, jsc_bnd, jec_bnd
call cpu_clock_begin(id_clock_forcing)
isc_bnd = index_bounds(1) ; iec_bnd = index_bounds(2)
jsc_bnd = index_bounds(3) ; jec_bnd = index_bounds(4)
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
Isq = G%IscB ; Ieq = G%IecB ; Jsq = G%JscB ; Jeq = G%JecB
isd = G%isd ; ied = G%ied ; jsd = G%jsd ; jed = G%jed
IsdB = G%IsdB ; IedB = G%IedB ; JsdB = G%JsdB ; JedB = G%JedB
isr = is-isd+1 ; ier = ie-isd+1 ; jsr = js-jsd+1 ; jer = je-jsd+1
i0 = is - isc_bnd ; j0 = js - jsc_bnd
kg_m2_s_conversion = US%kg_m2s_to_RZ_T
! allocation and initialization if this is the first time that this
! mechanical forcing type has been used.
if (.not.forces%initialized) then
call allocate_mech_forcing(G, forces, stress=.true., ustar=.not.CS%nonBous, &
press=.true., tau_mag=CS%nonBous)
call safe_alloc_ptr(forces%p_surf,isd,ied,jsd,jed)
call safe_alloc_ptr(forces%p_surf_full,isd,ied,jsd,jed)
if (CS%use_limited_P_SSH) then
forces%p_surf_SSH => forces%p_surf
else
forces%p_surf_SSH => forces%p_surf_full
endif
if (CS%rigid_sea_ice) then
call safe_alloc_ptr(forces%rigidity_ice_u,IsdB,IedB,jsd,jed)
call safe_alloc_ptr(forces%rigidity_ice_v,isd,ied,JsdB,JedB)
endif
forces%initialized = .true.
endif
if ( (associated(IOB%area_berg) .and. (.not. associated(forces%area_berg))) .or. &
(associated(IOB%mass_berg) .and. (.not. associated(forces%mass_berg))) ) &
call allocate_mech_forcing(G, forces, iceberg=.true.)
if (associated(IOB%ice_rigidity)) then
rigidity_at_h(:,:) = 0.0
call safe_alloc_ptr(forces%rigidity_ice_u,IsdB,IedB,jsd,jed)
call safe_alloc_ptr(forces%rigidity_ice_v,isd,ied,JsdB,JedB)
endif
forces%accumulate_rigidity = .true. ! Multiple components may contribute to rigidity.
if (associated(forces%rigidity_ice_u)) forces%rigidity_ice_u(:,:) = 0.0
if (associated(forces%rigidity_ice_v)) forces%rigidity_ice_v(:,:) = 0.0
! Set the weights for forcing fields that use running time averages.
if (present(reset_avg)) then ; if (reset_avg) forces%dt_force_accum = 0.0 ; endif
wt1 = 0.0 ; wt2 = 1.0
if (present(dt_forcing)) then
if ((forces%dt_force_accum > 0.0) .and. (dt_forcing > 0.0)) then
wt1 = forces%dt_force_accum / (forces%dt_force_accum + dt_forcing)
wt2 = 1.0 - wt1
endif
if (dt_forcing > 0.0) then
forces%dt_force_accum = max(forces%dt_force_accum, 0.0) + dt_forcing
else
forces%dt_force_accum = 0.0 ! Reset the averaging time interval.
endif
else
forces%dt_force_accum = 0.0 ! Reset the averaging time interval.
endif
! applied surface pressure from atmosphere and cryosphere
if (associated(IOB%p)) then
if (CS%max_p_surf >= 0.0) then
do j=js,je ; do i=is,ie
forces%p_surf_full(i,j) = G%mask2dT(i,j) * US%Pa_to_RL2_T2*IOB%p(i-i0,j-j0)
forces%p_surf(i,j) = MIN(forces%p_surf_full(i,j),CS%max_p_surf)
enddo ; enddo
else
do j=js,je ; do i=is,ie
forces%p_surf_full(i,j) = G%mask2dT(i,j) * US%Pa_to_RL2_T2*IOB%p(i-i0,j-j0)
forces%p_surf(i,j) = forces%p_surf_full(i,j)
enddo ; enddo
endif
else
do j=js,je ; do i=is,ie
forces%p_surf_full(i,j) = 0.0
forces%p_surf(i,j) = 0.0
enddo ; enddo
endif
forces%accumulate_p_surf = .true. ! Multiple components may contribute to surface pressure.
! Set the wind stresses and ustar.
if (wt1 <= 0.0) then
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, taux=forces%taux, tauy=forces%tauy, &
tau_halo=1)
if (associated(forces%ustar)) &
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, ustar=forces%ustar)
if (associated(forces%tau_mag)) &
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, mag_tau=forces%tau_mag)
else
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, taux=forces%taux, tauy=forces%tauy, &
tau_halo=1)
if (associated(forces%ustar)) then
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, ustar=ustar_tmp)
do j=js,je ; do i=is,ie
forces%ustar(i,j) = wt1*forces%ustar(i,j) + wt2*ustar_tmp(i,j)
enddo ; enddo
endif
if (associated(forces%tau_mag)) then
call extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, mag_tau=tau_mag_tmp)
do j=js,je ; do i=is,ie
forces%tau_mag(i,j) = wt1*forces%tau_mag(i,j) + wt2*tau_mag_tmp(i,j)
enddo ; enddo
endif
endif
! Find the net mass source in the input forcing without other adjustments.
if (CS%approx_net_mass_src .and. associated(forces%net_mass_src)) then
net_mass_src(:,:) = 0.0
i0 = is - isc_bnd ; j0 = js - jsc_bnd
do j=js,je ; do i=is,ie ; if (G%mask2dT(i,j) > 0.0) then
if (associated(IOB%lprec)) &
net_mass_src(i,j) = net_mass_src(i,j) + kg_m2_s_conversion * IOB%lprec(i-i0,j-j0)
if (associated(IOB%fprec)) &
net_mass_src(i,j) = net_mass_src(i,j) + kg_m2_s_conversion * IOB%fprec(i-i0,j-j0)
if (associated(IOB%seaice_melt)) &
net_mass_src(i,j) = net_mass_src(i,j) + kg_m2_s_conversion * IOB%seaice_melt(i-i0,j-j0)
if (associated(IOB%runoff)) &
net_mass_src(i,j) = net_mass_src(i,j) + kg_m2_s_conversion * IOB%runoff(i-i0,j-j0)
if (associated(IOB%calving)) &
net_mass_src(i,j) = net_mass_src(i,j) + kg_m2_s_conversion * IOB%calving(i-i0,j-j0)
if (associated(IOB%q_flux)) &
net_mass_src(i,j) = net_mass_src(i,j) - kg_m2_s_conversion * IOB%q_flux(i-i0,j-j0)
endif ; enddo ; enddo
if (wt1 <= 0.0) then
do j=js,je ; do i=is,ie
forces%net_mass_src(i,j) = wt2*net_mass_src(i,j)
enddo ; enddo
else
do j=js,je ; do i=is,ie
forces%net_mass_src(i,j) = wt1*forces%net_mass_src(i,j) + wt2*net_mass_src(i,j)
enddo ; enddo
endif
forces%net_mass_src_set = .true.
else
forces%net_mass_src_set = .false.
endif
! Obtain optional ice-berg related fluxes from the IOB type:
if (associated(IOB%area_berg)) then ; do j=js,je ; do i=is,ie
forces%area_berg(i,j) = IOB%area_berg(i-i0,j-j0) * G%mask2dT(i,j)
enddo ; enddo ; endif
if (associated(IOB%mass_berg)) then ; do j=js,je ; do i=is,ie
forces%mass_berg(i,j) = US%m_to_Z*US%kg_m3_to_R * IOB%mass_berg(i-i0,j-j0) * G%mask2dT(i,j)
enddo ; enddo ; endif
! Obtain sea ice related dynamic fields
if (associated(IOB%ice_rigidity)) then
do j=js,je ; do i=is,ie
rigidity_at_h(i,j) = US%m_to_L**3*US%Z_to_L*US%T_to_s * IOB%ice_rigidity(i-i0,j-j0) * G%mask2dT(i,j)
enddo ; enddo
call pass_var(rigidity_at_h, G%Domain, halo=1)
do I=is-1,ie ; do j=js,je
forces%rigidity_ice_u(I,j) = forces%rigidity_ice_u(I,j) + &
min(rigidity_at_h(i,j), rigidity_at_h(i+1,j))
enddo ; enddo
do i=is,ie ; do J=js-1,je
forces%rigidity_ice_v(i,J) = forces%rigidity_ice_v(i,J) + &
min(rigidity_at_h(i,j), rigidity_at_h(i,j+1))
enddo ; enddo
endif
if (CS%rigid_sea_ice) then
call pass_var(forces%p_surf_full, G%Domain, halo=1)
I_GEarth = 1.0 / CS%g_Earth
Kv_rho_ice = (CS%Kv_sea_ice / CS%density_sea_ice)
do I=is-1,ie ; do j=js,je
mass_ice = min(forces%p_surf_full(i,j), forces%p_surf_full(i+1,j)) * I_GEarth
mass_eff = 0.0
if (mass_ice > CS%rigid_sea_ice_mass) then
mass_eff = (mass_ice - CS%rigid_sea_ice_mass)**2 / (mass_ice + CS%rigid_sea_ice_mass)
endif
forces%rigidity_ice_u(I,j) = forces%rigidity_ice_u(I,j) + Kv_rho_ice * mass_eff
enddo ; enddo
do i=is,ie ; do J=js-1,je
mass_ice = min(forces%p_surf_full(i,j), forces%p_surf_full(i,j+1)) * I_GEarth
mass_eff = 0.0
if (mass_ice > CS%rigid_sea_ice_mass) then
mass_eff = (mass_ice - CS%rigid_sea_ice_mass)**2 / (mass_ice + CS%rigid_sea_ice_mass)
endif
forces%rigidity_ice_v(i,J) = forces%rigidity_ice_v(i,J) + Kv_rho_ice * mass_eff
enddo ; enddo
endif
if (CS%allow_flux_adjustments) then
! Apply adjustments to forces
call apply_force_adjustments(G, US, CS, Time, forces)
endif
!### ! Allow for user-written code to alter fluxes after all the above
!### call user_alter_mech_forcing(forces, Time, G, CS%urf_CS)
call cpu_clock_end(id_clock_forcing)
end subroutine convert_IOB_to_forces
!> This subroutine extracts the wind stresses and related fields like ustar from an
!! Ice_ocean_boundary_type into optional argument arrays, including changes of units, sign
!! conventions, and putting the fields into arrays with MOM-standard sized halos.
subroutine extract_IOB_stresses(IOB, index_bounds, Time, G, US, CS, taux, tauy, ustar, &
gustless_ustar, mag_tau, gustless_mag_tau, tau_halo)
type(ice_ocean_boundary_type), &
target, intent(in) :: IOB !< An ice-ocean boundary type with fluxes to drive
!! the ocean in a coupled model
integer, dimension(4), intent(in) :: index_bounds !< The i- and j- size of the arrays in IOB.
type(time_type), intent(in) :: Time !< The time of the fluxes, used for interpolating the
!! salinity to the right time, when it is being restored.
type(ocean_grid_type), intent(inout) :: G !< The ocean's grid structure
type(unit_scale_type), intent(in) :: US !< A dimensional unit scaling type
type(surface_forcing_CS),pointer :: CS !< A pointer to the control structure returned by a
!! previous call to surface_forcing_init.
real, dimension(SZIB_(G),SZJ_(G)), &
optional, intent(inout) :: taux !< The zonal wind stresses on a C-grid [R Z L T-2 ~> Pa].
real, dimension(SZI_(G),SZJB_(G)), &
optional, intent(inout) :: tauy !< The meridional wind stresses on a C-grid [R Z L T-2 ~> Pa].
real, dimension(SZI_(G),SZJ_(G)), &
optional, intent(inout) :: ustar !< The surface friction velocity [Z T-1 ~> m s-1].
real, dimension(SZI_(G),SZJ_(G)), &
optional, intent(out) :: gustless_ustar !< The surface friction velocity without
!! any contributions from gustiness [Z T-1 ~> m s-1].
real, dimension(SZI_(G),SZJ_(G)), &
optional, intent(inout) :: mag_tau !< The magintude of the wind stress at tracer points
!! including subgridscale variability and gustiness [R Z L T-2 ~> Pa]
real, dimension(SZI_(G),SZJ_(G)), &
optional, intent(out) :: gustless_mag_tau !< The magintude of the wind stress at tracer points
!! without any contributions from gustiness [R Z L T-2 ~> Pa]
integer, optional, intent(in) :: tau_halo !< The halo size of wind stresses to set, 0 by default.
! Local variables
real, dimension(SZI_(G),SZJ_(G)) :: taux_in_A ! Zonal wind stresses [R Z L T-2 ~> Pa] at h points
real, dimension(SZI_(G),SZJ_(G)) :: tauy_in_A ! Meridional wind stresses [R Z L T-2 ~> Pa] at h points
real, dimension(SZIB_(G),SZJ_(G)) :: taux_in_C ! Zonal wind stresses [R Z L T-2 ~> Pa] at u points
real, dimension(SZI_(G),SZJB_(G)) :: tauy_in_C ! Meridional wind stresses [R Z L T-2 ~> Pa] at v points
real, dimension(SZIB_(G),SZJB_(G)) :: taux_in_B ! Zonal wind stresses [R Z L T-2 ~> Pa] at q points
real, dimension(SZIB_(G),SZJB_(G)) :: tauy_in_B ! Meridional wind stresses [R Z L T-2 ~> Pa] at q points
real :: gustiness ! unresolved gustiness that contributes to ustar [R Z L T-2 ~> Pa]
real :: Irho0 ! Inverse of the mean density rescaled to [Z L-1 R-1 ~> m3 kg-1]
real :: taux2, tauy2 ! squared wind stresses [R2 Z2 L2 T-4 ~> Pa2]
real :: tau_mag ! magnitude of the wind stress [R Z L T-2 ~> Pa]
real :: stress_conversion ! A unit conversion factor from Pa times any stress multiplier [R Z L T-2 Pa-1 ~> 1]
logical :: do_ustar, do_gustless, do_tau_mag, do_gustless_tau_mag
integer :: wind_stagger ! AGRID, BGRID_NE, or CGRID_NE (integers from MOM_domains)
integer :: i, j, is, ie, js, je, ish, ieh, jsh, jeh, Isqh, Ieqh, Jsqh, Jeqh, i0, j0, halo
halo = 0 ; if (present(tau_halo)) halo = tau_halo
is = G%isc ; ie = G%iec ; js = G%jsc ; je = G%jec
ish = G%isc-halo ; ieh = G%iec+halo ; jsh = G%jsc-halo ; jeh = G%jec+halo
Isqh = G%IscB-halo ; Ieqh = G%IecB+halo ; Jsqh = G%JscB-halo ; Jeqh = G%JecB+halo
i0 = is - index_bounds(1) ; j0 = js - index_bounds(3)
IRho0 = US%L_to_Z / CS%Rho0
stress_conversion = US%Pa_to_RLZ_T2 * CS%wind_stress_multiplier
do_ustar = present(ustar) ; do_gustless = present(gustless_ustar)
do_tau_mag = present(mag_tau) ; do_gustless_tau_mag = present(gustless_mag_tau)
wind_stagger = CS%wind_stagger
if ((IOB%wind_stagger == AGRID) .or. (IOB%wind_stagger == BGRID_NE) .or. &
(IOB%wind_stagger == CGRID_NE)) wind_stagger = IOB%wind_stagger
if (associated(IOB%u_flux).neqv.associated(IOB%v_flux)) call MOM_error(FATAL,"extract_IOB_stresses: "//&
"associated(IOB%u_flux) /= associated(IOB%v_flux !!!")
if (present(taux).neqv.present(tauy)) call MOM_error(FATAL,"extract_IOB_stresses: "//&
"present(taux) /= present(tauy) !!!")
! Set surface momentum stress related fields as a function of staggering.
if (present(taux) .or. present(tauy) .or. &
((do_ustar .or. do_tau_mag .or. do_gustless .or. do_gustless_tau_mag) &
.and. .not.associated(IOB%stress_mag)) ) then
if (wind_stagger == BGRID_NE) then
taux_in_B(:,:) = 0.0 ; tauy_in_B(:,:) = 0.0
if (associated(IOB%u_flux).and.associated(IOB%v_flux)) then
do J=js,je ; do I=is,ie
taux_in_B(I,J) = IOB%u_flux(i-i0,j-j0) * stress_conversion