forked from mom-ocean/MOM6
-
Notifications
You must be signed in to change notification settings - Fork 65
/
Copy pathMOM_neutral_diffusion.F90
2576 lines (2326 loc) · 137 KB
/
MOM_neutral_diffusion.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
!> A column-wise toolbox for implementing neutral diffusion
module MOM_neutral_diffusion
! This file is part of MOM6. See LICENSE.md for the license.
use MOM_cpu_clock, only : cpu_clock_id, cpu_clock_begin, cpu_clock_end
use MOM_cpu_clock, only : CLOCK_MODULE, CLOCK_ROUTINE
use MOM_diag_mediator, only : diag_ctrl, time_type
use MOM_diag_mediator, only : post_data, register_diag_field
use MOM_EOS, only : EOS_type, EOS_manual_init, calculate_compress, calculate_density_derivs
use MOM_EOS, only : calculate_density_second_derivs
use MOM_EOS, only : extract_member_EOS, EOS_LINEAR, EOS_TEOS10, EOS_WRIGHT
use MOM_error_handler, only : MOM_error, FATAL, WARNING, MOM_mesg, is_root_pe
use MOM_file_parser, only : get_param, log_version, param_file_type
use MOM_file_parser, only : openParameterBlock, closeParameterBlock
use MOM_grid, only : ocean_grid_type
use MOM_remapping, only : remapping_CS, initialize_remapping
use MOM_remapping, only : extract_member_remapping_CS, build_reconstructions_1d
use MOM_remapping, only : average_value_ppoly, remappingSchemesDoc, remappingDefaultScheme
use MOM_tracer_registry, only : tracer_registry_type
use MOM_verticalGrid, only : verticalGrid_type
use polynomial_functions, only : evaluation_polynomial, first_derivative_polynomial
use PPM_functions, only : PPM_reconstruction, PPM_boundary_extrapolation
use regrid_edge_values, only : edge_values_implicit_h4
implicit none ; private
#include <MOM_memory.h>
public neutral_diffusion
public neutral_diffusion_init
public neutral_diffusion_diag_init
public neutral_diffusion_end
public neutral_diffusion_calc_coeffs
public neutral_diffusion_unit_tests
type, public :: neutral_diffusion_CS ; private
integer :: nkp1 ! Number of interfaces for a column = nk + 1
integer :: nsurf ! Number of neutral surfaces
integer :: ppoly_deg = 2 ! Degree of polynomial used for reconstructions
logical :: continuous_reconstruction = .true. ! True if using continuous PPM reconstruction at interfaces
logical :: boundary_extrap = .true.
logical :: refine_position = .false.
integer :: max_iter ! Maximum number of iterations if refine_position is defined
real :: tolerance ! Convergence criterion representing difference from true neutrality
real :: ref_pres ! Reference pressure, negative if using locally referenced neutral density
! Positions of neutral surfaces in both the u, v directions
real, allocatable, dimension(:,:,:) :: uPoL ! Non-dimensional position with left layer uKoL-1, u-point
real, allocatable, dimension(:,:,:) :: uPoR ! Non-dimensional position with right layer uKoR-1, u-point
integer, allocatable, dimension(:,:,:) :: uKoL ! Index of left interface corresponding to neutral surface, u-point
integer, allocatable, dimension(:,:,:) :: uKoR ! Index of right interface corresponding to neutral surface, u-point
real, allocatable, dimension(:,:,:) :: uHeff ! Effective thickness at u-point (H units)
real, allocatable, dimension(:,:,:) :: vPoL ! Non-dimensional position with left layer uKoL-1, v-point
real, allocatable, dimension(:,:,:) :: vPoR ! Non-dimensional position with right layer uKoR-1, v-point
integer, allocatable, dimension(:,:,:) :: vKoL ! Index of left interface corresponding to neutral surface, v-point
integer, allocatable, dimension(:,:,:) :: vKoR ! Index of right interface corresponding to neutral surface, v-point
real, allocatable, dimension(:,:,:) :: vHeff ! Effective thickness at v-point (H units)
! Coefficients of polynomial reconstructions for temperature and salinity
real, allocatable, dimension(:,:,:,:) :: ppoly_coeffs_T !< Polynomial coefficients for temperature
real, allocatable, dimension(:,:,:,:) :: ppoly_coeffs_S !< Polynomial coefficients for temperature
! Variables needed for continuous reconstructions
real, allocatable, dimension(:,:,:) :: dRdT ! dRho/dT (kg/m3/degC) at interfaces
real, allocatable, dimension(:,:,:) :: dRdS ! dRho/dS (kg/m3/ppt) at interfaces
real, allocatable, dimension(:,:,:) :: Tint ! Interface T (degC)
real, allocatable, dimension(:,:,:) :: Sint ! Interface S (ppt)
real, allocatable, dimension(:,:,:) :: Pint ! Interface pressure (Pa)
! Variables needed for discontinuous reconstructions
real, allocatable, dimension(:,:,:,:) :: T_i ! Top edge reconstruction of temperature (degC)
real, allocatable, dimension(:,:,:,:) :: S_i ! Top edge reconstruction of salinity (ppt)
real, allocatable, dimension(:,:,:,:) :: dRdT_i ! dRho/dT (kg/m3/degC) at top edge
real, allocatable, dimension(:,:,:,:) :: dRdS_i ! dRho/dS (kg/m3/ppt) at top edge
type(diag_ctrl), pointer :: diag ! structure to regulate output
integer, allocatable, dimension(:) :: id_neutral_diff_tracer_conc_tend ! tracer concentration tendency
integer, allocatable, dimension(:) :: id_neutral_diff_tracer_cont_tend ! tracer content tendency
integer, allocatable, dimension(:) :: id_neutral_diff_tracer_cont_tend_2d ! k-summed tracer content tendency
integer, allocatable, dimension(:) :: id_neutral_diff_tracer_trans_x_2d ! k-summed ndiff zonal tracer transport
integer, allocatable, dimension(:) :: id_neutral_diff_tracer_trans_y_2d ! k-summed ndiff merid tracer transport
real :: C_p ! heat capacity of seawater (J kg-1 K-1)
type(remapping_CS) :: remap_CS
end type neutral_diffusion_CS
! This include declares and sets the variable "version".
#include "version_variable.h"
character(len=40) :: mdl = "MOM_neutral_diffusion" ! module name
logical :: debug_this_module = .false. ! If true, verbose output of find neutral position
contains
!> Read parameters and allocate control structure for neutral_diffusion module.
logical function neutral_diffusion_init(Time, G, param_file, diag, CS)
type(time_type), target, intent(in) :: Time !< Time structure
type(ocean_grid_type), intent(in) :: G !< Grid structure
type(diag_ctrl), target, intent(inout) :: diag !< Diagnostics control structure
type(param_file_type), intent(in) :: param_file !< Parameter file structure
type(neutral_diffusion_CS), pointer :: CS !< Neutral diffusion control structure
! Local variables
character(len=256) :: mesg ! Message for error messages.
character(len=80) :: string ! Temporary strings
if (associated(CS)) then
call MOM_error(FATAL, "neutral_diffusion_init called with associated control structure.")
return
endif
! Log this module and master switch for turning it on/off
call log_version(param_file, mdl, version, &
"This module implements neutral diffusion of tracers")
call get_param(param_file, mdl, "USE_NEUTRAL_DIFFUSION", neutral_diffusion_init, &
"If true, enables the neutral diffusion module.", &
default=.false.)
if (.not.neutral_diffusion_init) then
return
endif
allocate(CS)
CS%diag => diag
! call openParameterBlock(param_file,'NEUTRAL_DIFF')
! Read all relevant parameters and write them to the model log.
call get_param(param_file, mdl, "NDIFF_CONTINUOUS", CS%continuous_reconstruction, &
"If true, uses a continuous reconstruction of T and S when \n"// &
"finding neutral surfaces along which diffusion will happen.\n"// &
"If false, a PPM discontinuous reconstruction of T and S \n"// &
"is done which results in a higher order routine but exacts \n"// &
"a higher computational cost.", default=.true.)
call get_param(param_file, mdl, "NDIFF_REF_PRES", CS%ref_pres, &
"The reference pressure (Pa) used for the derivatives of \n"// &
"the equation of state. If negative (default), local \n"// &
"pressure is used.", &
default = -1.)
! Initialize and configure remapping
if (CS%continuous_reconstruction .eqv. .false.) then
call get_param(param_file, mdl, "NDIFF_BOUNDARY_EXTRAP", CS%boundary_extrap, &
"Uses a rootfinding approach to find the position of a\n"// &
"neutral surface within a layer taking into account the\n"// &
"nonlinearity of the equation of state and the\n"// &
"polynomial reconstructions of T/S.", &
default=.false.)
call get_param(param_file, mdl, "NDIFF_REMAPPING_SCHEME", string, &
"This sets the reconstruction scheme used\n"//&
"for vertical remapping for all variables.\n"//&
"It can be one of the following schemes:\n"//&
trim(remappingSchemesDoc), default=remappingDefaultScheme)
call initialize_remapping( CS%remap_CS, string, boundary_extrapolation = CS%boundary_extrap )
call extract_member_remapping_CS(CS%remap_CS, degree=CS%ppoly_deg)
call get_param(param_file, mdl, "NDIFF_REFINE_POSITION", CS%refine_position, &
"Uses a rootfinding approach to find the position of a\n"// &
"neutral surface within a layer taking into account the\n"// &
"nonlinearity of the equation of state and the\n"// &
"polynomial reconstructions of T/S.", &
default=.false.)
if (CS%refine_position) then
call get_param(param_file, mdl, "NDIFF_TOLERANCE", CS%tolerance, &
"Sets the convergence criterion for finding the neutral\n"// &
"position within a layer in kg m-3.", &
default=1.e-10)
call get_param(param_file, mdl, "NDIFF_MAX_ITER", CS%max_iter, &
"The maximum number of iterations to be done before \n"// &
"exiting the iterative loop to find the neutral surface", &
default=10)
endif
call get_param(param_file, mdl, "NDIFF_DEBUG", debug_this_module, &
"Turns on verbose output for discontinuous neutral \n"// &
"diffusion routines.", &
default = .false.)
endif
! call get_param(param_file, mdl, "KHTR", CS%KhTr, &
! "The background along-isopycnal tracer diffusivity.", &
! units="m2 s-1", default=0.0)
! call closeParameterBlock(param_file)
if (CS%continuous_reconstruction) then
CS%nsurf = 2*G%ke+2 ! Continuous reconstruction means that every interface has two connections
allocate(CS%dRdT(SZI_(G),SZJ_(G),SZK_(G)+1)) ; CS%dRdT(:,:,:) = 0.
allocate(CS%dRdS(SZI_(G),SZJ_(G),SZK_(G)+1)) ; CS%dRdS(:,:,:) = 0.
else
CS%nsurf = 4*G%ke ! Discontinuous means that every interface has four connections
allocate(CS%T_i(SZI_(G),SZJ_(G),SZK_(G),2)) ; CS%T_i(:,:,:,:) = 0.
allocate(CS%S_i(SZI_(G),SZJ_(G),SZK_(G),2)) ; CS%S_i(:,:,:,:) = 0.
allocate(CS%dRdT_i(SZI_(G),SZJ_(G),SZK_(G),2)) ; CS%dRdT_i(:,:,:,:) = 0.
allocate(CS%dRdS_i(SZI_(G),SZJ_(G),SZK_(G),2)) ; CS%dRdS_i(:,:,:,:) = 0.
allocate(CS%ppoly_coeffs_T(SZI_(G),SZJ_(G),SZK_(G),CS%ppoly_deg+1)) ; CS%ppoly_coeffs_T(:,:,:,:) = 0.
allocate(CS%ppoly_coeffs_S(SZI_(G),SZJ_(G),SZK_(G),CS%ppoly_deg+1)) ; CS%ppoly_coeffs_S(:,:,:,:) = 0.
endif
! T-points
allocate(CS%Tint(SZI_(G),SZJ_(G),SZK_(G)+1)) ; CS%Tint(:,:,:) = 0.
allocate(CS%Sint(SZI_(G),SZJ_(G),SZK_(G)+1)) ; CS%Sint(:,:,:) = 0.
allocate(CS%Pint(SZI_(G),SZJ_(G),SZK_(G)+1)) ; CS%Pint(:,:,:) = 0.
! U-points
allocate(CS%uPoL(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%uPoL(G%isc-1:G%iec,G%jsc:G%jec,:) = 0.
allocate(CS%uPoR(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%uPoR(G%isc-1:G%iec,G%jsc:G%jec,:) = 0.
allocate(CS%uKoL(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%uKoL(G%isc-1:G%iec,G%jsc:G%jec,:) = 0
allocate(CS%uKoR(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%uKoR(G%isc-1:G%iec,G%jsc:G%jec,:) = 0
allocate(CS%uHeff(G%isd:G%ied,G%jsd:G%jed,CS%nsurf-1)); CS%uHeff(G%isc-1:G%iec,G%jsc:G%jec,:) = 0
! V-points
allocate(CS%vPoL(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%vPoL(G%isc:G%iec,G%jsc-1:G%jec,:) = 0.
allocate(CS%vPoR(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%vPoR(G%isc:G%iec,G%jsc-1:G%jec,:) = 0.
allocate(CS%vKoL(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%vKoL(G%isc:G%iec,G%jsc-1:G%jec,:) = 0
allocate(CS%vKoR(G%isd:G%ied,G%jsd:G%jed, CS%nsurf)); CS%vKoR(G%isc:G%iec,G%jsc-1:G%jec,:) = 0
allocate(CS%vHeff(G%isd:G%ied,G%jsd:G%jed,CS%nsurf-1)); CS%vHeff(G%isc:G%iec,G%jsc-1:G%jec,:) = 0
end function neutral_diffusion_init
!> Diagnostic handles for neutral diffusion tendencies.
subroutine neutral_diffusion_diag_init(Time, G, diag, C_p, Reg, CS)
type(time_type),target, intent(in) :: Time !< Time structure
type(ocean_grid_type), intent(in) :: G !< Grid structure
type(diag_ctrl), intent(in) :: diag !< Diagnostics control structure
type(tracer_registry_type), intent(in) :: Reg !< Tracer structure
real, intent(in) :: C_p !< Seawater heat capacity
type(neutral_diffusion_CS), pointer :: CS !< Neutral diffusion control structure
! local
integer :: n,ntr
if(.not. associated(CS)) return
ntr = Reg%ntr
CS%C_p = C_p
allocate(CS%id_neutral_diff_tracer_conc_tend(ntr))
allocate(CS%id_neutral_diff_tracer_cont_tend(ntr))
allocate(CS%id_neutral_diff_tracer_cont_tend_2d(ntr))
allocate(CS%id_neutral_diff_tracer_trans_x_2d(ntr))
allocate(CS%id_neutral_diff_tracer_trans_y_2d(ntr))
CS%id_neutral_diff_tracer_conc_tend(:) = -1
CS%id_neutral_diff_tracer_cont_tend(:) = -1
CS%id_neutral_diff_tracer_cont_tend_2d(:) = -1
CS%id_neutral_diff_tracer_trans_x_2d(:) = -1
CS%id_neutral_diff_tracer_trans_y_2d(:) = -1
do n=1,ntr
if(trim(Reg%Tr(n)%name) == 'T') then
CS%id_neutral_diff_tracer_conc_tend(n) = register_diag_field('ocean_model', &
'ndiff_tracer_conc_tendency_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Neutral diffusion tracer concentration tendency for '//trim(Reg%Tr(n)%name),&
'Degree C per second')
CS%id_neutral_diff_tracer_cont_tend(n) = register_diag_field('ocean_model', &
'ndiff_tracer_cont_tendency_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Neutral diffusion tracer content tendency for '//trim(Reg%Tr(n)%name), &
'Watts/m2',cmor_field_name='opottemppmdiff', cmor_units='W m-2', &
cmor_standard_name= &
'tendency_of_sea_water_potential_temperature_expressed_as_heat_content_due_to_parameterized_mesocale_diffusion', &
cmor_long_name = &
'Tendency of sea water potential temperature expressed as heat content due to parameterized mesocale diffusion', &
v_extensive=.true.)
CS%id_neutral_diff_tracer_cont_tend_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_cont_tendency_2d_'//trim(Reg%Tr(n)%name), diag%axesT1, Time, &
'Depth integrated neutral diffusion tracer content tendency for '//trim(Reg%Tr(n)%name), &
'Watts/m2',cmor_field_name='opottemppmdiff_2d', cmor_units='W m-2', &
cmor_standard_name= &
'tendency_of_sea_water_potential_temperature_expressed_as_heat_content_due_to_parameterized_mesocale_diffusion_depth_integrated',&
cmor_long_name = &
'Tendency of sea water potential temperature expressed as heat content due to parameterized mesocale diffusion depth integrated')
CS%id_neutral_diff_tracer_trans_x_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_trans_x_2d_'//trim(Reg%Tr(n)%name), diag%axesCu1, Time, &
'Depth integrated neutral diffusion zonal tracer transport for '//trim(Reg%Tr(n)%name),&
'Watts')
CS%id_neutral_diff_tracer_trans_y_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_trans_y_2d_'//trim(Reg%Tr(n)%name), diag%axesCv1, Time, &
'Depth integrated neutral diffusion merid tracer transport for '//trim(Reg%Tr(n)%name),&
'Watts')
elseif(trim(Reg%Tr(n)%name) == 'S') then
CS%id_neutral_diff_tracer_conc_tend(n) = register_diag_field('ocean_model', &
'ndiff_tracer_conc_tendency_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Neutral diffusion tracer concentration tendency for '//trim(Reg%Tr(n)%name),&
'tracer concentration units per second')
CS%id_neutral_diff_tracer_cont_tend(n) = register_diag_field('ocean_model', &
'ndiff_tracer_cont_tendency_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Neutral diffusion tracer content tendency for '//trim(Reg%Tr(n)%name), &
'kg m-2 s-1',cmor_field_name='osaltpmdiff', cmor_units='kg m-2 s-1', &
cmor_standard_name= &
'tendency_of_sea_water_salinity_expressed_as_salt_content_due_to_parameterized_mesocale_diffusion', &
cmor_long_name = &
'Tendency of sea water salinity expressed as salt content due to parameterized mesocale diffusion', &
v_extensive=.true.)
CS%id_neutral_diff_tracer_cont_tend_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_cont_tendency_2d_'//trim(Reg%Tr(n)%name), diag%axesT1, Time, &
'Depth integrated neutral diffusion tracer content tendency for '//trim(Reg%Tr(n)%name), &
'kg m-2 s-1',cmor_field_name='osaltpmdiff_2d', cmor_units='kg m-2 s-1', &
cmor_standard_name= &
'tendency_of_sea_water_salinity_expressed_as_salt_content_due_to_parameterized_mesocale_diffusion_depth_integrated',&
cmor_long_name = &
'Tendency of sea water salinity expressed as salt content due to parameterized mesocale diffusion depth integrated')
CS%id_neutral_diff_tracer_trans_x_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_trans_x_2d_'//trim(Reg%Tr(n)%name), diag%axesCu1, Time, &
'Depth integrated neutral diffusion zonal tracer transport for '//trim(Reg%Tr(n)%name),&
'kg/s')
CS%id_neutral_diff_tracer_trans_y_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_trans_y_2d_'//trim(Reg%Tr(n)%name), diag%axesCv1, Time, &
'Depth integrated neutral diffusion merid tracer transport for '//trim(Reg%Tr(n)%name),&
'kg/s')
else
CS%id_neutral_diff_tracer_conc_tend(n) = register_diag_field('ocean_model', &
'ndiff_tracer_conc_tendency_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Neutral diffusion tracer concentration tendency for '//trim(Reg%Tr(n)%name),&
'tracer concentration * m-2 s-1')
CS%id_neutral_diff_tracer_cont_tend(n) = register_diag_field('ocean_model',&
'ndiff_tracer_cont_tendency_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Neutral diffusion tracer content tendency for '//trim(Reg%Tr(n)%name), &
'tracer content * m-2 s-1', v_extensive=.true.)
CS%id_neutral_diff_tracer_cont_tend_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_cont_tendency_2d_'//trim(Reg%Tr(n)%name), diag%axesTL, Time, &
'Depth integrated neutral diffusion tracer content tendency for '//trim(Reg%Tr(n)%name),&
'tracer content * m-2 s-1')
CS%id_neutral_diff_tracer_trans_x_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_trans_x_2d_'//trim(Reg%Tr(n)%name), diag%axesCu1, Time, &
'Depth integrated neutral diffusion zonal tracer transport for '//trim(Reg%Tr(n)%name),&
'kg/s')
CS%id_neutral_diff_tracer_trans_y_2d(n) = register_diag_field('ocean_model', &
'ndiff_tracer_trans_y_2d_'//trim(Reg%Tr(n)%name), diag%axesCv1, Time, &
'Depth integrated neutral diffusion merid tracer transport for '//trim(Reg%Tr(n)%name),&
'kg/s')
endif
enddo
end subroutine neutral_diffusion_diag_init
!> Calculate remapping factors for u/v columns used to map adjoining columns to
!! a shared coordinate space.
subroutine neutral_diffusion_calc_coeffs(G, GV, h, T, S, EOS, CS)
type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness (H units)
real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: T !< Potential temperature (degC)
real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: S !< Salinity (ppt)
type(EOS_type), pointer :: EOS !< Equation of state structure
type(neutral_diffusion_CS), pointer :: CS !< Neutral diffusion control structure
! Local variables
integer :: i, j, k
! Variables used for reconstructions
real, dimension(SZK_(G),2) :: ppoly_r_S ! Reconstruction slopes
integer :: iMethod
real, dimension(SZI_(G)+1) :: ref_pres ! Reference pressure used to calculate alpha/beta
! If doing along isopycnal diffusion (as opposed to neutral diffusion, set the reference pressure)
if (CS%ref_pres>=0.) ref_pres(:) = CS%ref_pres
if (CS%continuous_reconstruction) then
CS%dRdT(:,:,:) = 0.
CS%dRdS(:,:,:) = 0.
else
CS%T_i(:,:,:,:) = 0.
CS%S_i(:,:,:,:) = 0.
CS%dRdT_i(:,:,:,:) = 0.
CS%dRdS_i(:,:,:,:) = 0.
endif
! Calculate pressure at interfaces
CS%Pint(:,:,1) = 0.
do k=1,G%ke ; do j=G%jsc-1, G%jec+1 ; do i=G%isc-1,G%iec+1
CS%Pint(i,j,k+1) = CS%Pint(i,j,k) + h(i,j,k)*GV%H_to_Pa
enddo ; enddo ; enddo
do j = G%jsc-1, G%jec+1
! Interpolate state to interface
do i = G%isc-1, G%iec+1
if (CS%continuous_reconstruction) then
call interface_scalar(G%ke, h(i,j,:), T(i,j,:), CS%Tint(i,j,:), 2)
call interface_scalar(G%ke, h(i,j,:), S(i,j,:), CS%Sint(i,j,:), 2)
else
call build_reconstructions_1d( CS%remap_CS, G%ke, h(i,j,:), T(i,j,:), CS%ppoly_coeffs_T(i,j,:,:), &
CS%T_i(i,j,:,:), ppoly_r_S, iMethod )
call build_reconstructions_1d( CS%remap_CS, G%ke, h(i,j,:), S(i,j,:), CS%ppoly_coeffs_S(i,j,:,:), &
CS%S_i(i,j,:,:), ppoly_r_S, iMethod )
endif
enddo
! Continuous reconstruction
if (CS%continuous_reconstruction) then
do k = 1, G%ke+1
if (CS%ref_pres<0) ref_pres(:) = CS%Pint(:,j,k)
call calculate_density_derivs(CS%Tint(:,j,k), CS%Sint(:,j,k), ref_pres, &
CS%dRdT(:,j,k), CS%dRdS(:,j,k), G%isc-1, G%iec-G%isc+3, EOS)
enddo
else ! Discontinuous reconstruction
do k = 1, G%ke
if (CS%ref_pres<0) ref_pres(:) = CS%Pint(:,j,k)
call calculate_density_derivs(CS%T_i(:,j,k,1), CS%S_i(:,j,k,1), ref_pres, &
CS%dRdT_i(:,j,k,1), CS%dRdS_i(:,j,k,1), G%isc-1, G%iec-G%isc+3, EOS)
if (CS%ref_pres<0) ref_pres(:) = CS%Pint(:,j,k+1)
call calculate_density_derivs(CS%T_i(:,j,k,2), CS%S_i(:,j,k,2), ref_pres, &
CS%dRdT_i(:,j,k,2), CS%dRdS_i(:,j,k,2), G%isc-1, G%iec-G%isc+3, EOS)
enddo
endif
enddo
CS%uhEff(:,:,:) = 0.
CS%vhEff(:,:,:) = 0.
CS%uPoL(:,:,:) = 0.
CS%vPoL(:,:,:) = 0.
CS%uPoR(:,:,:) = 0.
CS%vPoR(:,:,:) = 0.
CS%uKoL(:,:,:) = 1
CS%vKoL(:,:,:) = 1
CS%uKoR(:,:,:) = 1
CS%vKoR(:,:,:) = 1
! Neutral surface factors at U points
do j = G%jsc, G%jec ; do I = G%isc-1, G%iec
if (G%mask2dCu(I,j) > 0.) then
if (CS%continuous_reconstruction) then
call find_neutral_surface_positions_continuous(G%ke, &
CS%Pint(i,j,:), CS%Tint(i,j,:), CS%Sint(i,j,:), CS%dRdT(i,j,:), CS%dRdS(i,j,:), &
CS%Pint(i+1,j,:), CS%Tint(i+1,j,:), CS%Sint(i+1,j,:), CS%dRdT(i+1,j,:), CS%dRdS(i+1,j,:), &
CS%uPoL(I,j,:), CS%uPoR(I,j,:), CS%uKoL(I,j,:), CS%uKoR(I,j,:), CS%uhEff(I,j,:) )
else
call find_neutral_surface_positions_discontinuous(G%ke, CS%ppoly_deg, &
CS%Pint(i,j,:), CS%T_i(i,j,:,:), CS%S_i(i,j,:,:), CS%dRdT_i(i,j,:,:), CS%dRdS_i(i,j,:,:), &
CS%Pint(i+1,j,:), CS%T_i(i+1,j,:,:), CS%S_i(i+1,j,:,:), CS%dRdT_i(i+1,j,:,:), CS%dRdS_i(i+1,j,:,:), &
CS%uPoL(I,j,:), CS%uPoR(I,j,:), CS%uKoL(I,j,:), CS%uKoR(I,j,:), CS%uhEff(I,j,:), &
CS%refine_position, CS%ppoly_coeffs_T(i,j,:,:), CS%ppoly_coeffs_S(i,j,:,:), &
CS%ppoly_coeffs_T(i+1,j,:,:), CS%ppoly_coeffs_S(i+1,j,:,:), EOS, CS%max_iter, CS%tolerance, CS%ref_pres)
endif
endif
enddo ; enddo
! Neutral surface factors at V points
do J = G%jsc-1, G%jec ; do i = G%isc, G%iec
if (G%mask2dCv(i,J) > 0.) then
if (CS%continuous_reconstruction) then
call find_neutral_surface_positions_continuous(G%ke, &
CS%Pint(i,j,:), CS%Tint(i,j,:), CS%Sint(i,j,:), CS%dRdT(i,j,:), CS%dRdS(i,j,:), &
CS%Pint(i,j+1,:), CS%Tint(i,j+1,:), CS%Sint(i,j+1,:), CS%dRdT(i,j+1,:), CS%dRdS(i,j+1,:), &
CS%vPoL(i,J,:), CS%vPoR(i,J,:), CS%vKoL(i,J,:), CS%vKoR(i,J,:), CS%vhEff(i,J,:) )
else
call find_neutral_surface_positions_discontinuous(G%ke, CS%ppoly_deg, &
CS%Pint(i,j,:), CS%T_i(i,j,:,:), CS%S_i(i,j,:,:), CS%dRdT_i(i,j,:,:), CS%dRdS_i(i,j,:,:), &
CS%Pint(i,j+1,:), CS%T_i(i,j+1,:,:), CS%S_i(i,j+1,:,:), CS%dRdT_i(i,j+1,:,:), CS%dRdS_i(i,j+1,:,:), &
CS%vPoL(I,j,:), CS%vPoR(I,j,:), CS%vKoL(I,j,:), CS%vKoR(I,j,:), CS%vhEff(I,j,:), &
CS%refine_position, CS%ppoly_coeffs_T(i,j,:,:), CS%ppoly_coeffs_S(i,j,:,:), &
CS%ppoly_coeffs_T(i,j+1,:,:), CS%ppoly_coeffs_S(i,j+1,:,:), EOS, CS%max_iter, CS%tolerance, CS%ref_pres)
endif
endif
enddo ; enddo
CS%uhEff(:,:,:) = CS%uhEff(:,:,:) / GV%H_to_pa
CS%vhEff(:,:,:) = CS%vhEff(:,:,:) / GV%H_to_pa
end subroutine neutral_diffusion_calc_coeffs
!> Update tracer concentration due to neutral diffusion; layer thickness unchanged by this update.
subroutine neutral_diffusion(G, GV, h, Coef_x, Coef_y, Tracer, m, dt, name, CS)
type(ocean_grid_type), intent(in) :: G !< Ocean grid structure
type(verticalGrid_type), intent(in) :: GV !< ocean vertical grid structure
real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(in) :: h !< Layer thickness (H units)
real, dimension(SZIB_(G),SZJ_(G)), intent(in) :: Coef_x !< dt * Kh * dy / dx at u-points (m^2)
real, dimension(SZI_(G),SZJB_(G)), intent(in) :: Coef_y !< dt * Kh * dx / dy at u-points (m^2)
real, dimension(SZI_(G),SZJ_(G),SZK_(G)), intent(inout) :: Tracer !< Tracer concentration
integer, intent(in) :: m !< Tracer number
real, intent(in) :: dt !< Tracer time step * I_numitts (I_numitts in tracer_hordiff)
character(len=32), intent(in) :: name !< Tracer name
type(neutral_diffusion_CS), pointer :: CS !< Neutral diffusion control structure
! Local variables
real, dimension(SZIB_(G),SZJ_(G),CS%nsurf-1) :: uFlx ! Zonal flux of tracer (concentration * H)
real, dimension(SZI_(G),SZJB_(G),CS%nsurf-1) :: vFlx ! Meridional flux of tracer (concentration * H)
real, dimension(SZI_(G),SZJ_(G),G%ke) :: tendency ! tendency array for diagn
real, dimension(SZI_(G),SZJ_(G)) :: tendency_2d ! depth integrated content tendency for diagn
real, dimension(SZIB_(G),SZJ_(G)) :: trans_x_2d ! depth integrated diffusive tracer x-transport diagn
real, dimension(SZI_(G),SZJB_(G)) :: trans_y_2d ! depth integrated diffusive tracer y-transport diagn
real, dimension(G%ke) :: dTracer ! change in tracer concentration due to ndiffusion
integer :: i, j, k, ks, nk
real :: ppt2mks, Idt, convert
nk = GV%ke
! for diagnostics
if(CS%id_neutral_diff_tracer_conc_tend(m) > 0 .or. &
CS%id_neutral_diff_tracer_cont_tend(m) > 0 .or. &
CS%id_neutral_diff_tracer_cont_tend_2d(m) > 0 .or. &
CS%id_neutral_diff_tracer_trans_x_2d(m) > 0 .or. &
CS%id_neutral_diff_tracer_trans_y_2d(m) > 0) then
ppt2mks = 0.001
Idt = 1.0/dt
tendency(:,:,:) = 0.0
tendency_2d(:,:) = 0.0
trans_x_2d(:,:) = 0.0
trans_y_2d(:,:) = 0.0
convert = 1.0
if(trim(name) == 'T') convert = CS%C_p * GV%H_to_kg_m2
if(trim(name) == 'S') convert = ppt2mks * GV%H_to_kg_m2
endif
uFlx(:,:,:) = 0.
vFlx(:,:,:) = 0.
! x-flux
do j = G%jsc,G%jec ; do I = G%isc-1,G%iec
if (G%mask2dCu(I,j)>0.) then
call neutral_surface_flux(nk, CS%nsurf, CS%ppoly_deg, h(i,j,:), h(i+1,j,:), &
Tracer(i,j,:), Tracer(i+1,j,:), &
CS%uPoL(I,j,:), CS%uPoR(I,j,:), &
CS%uKoL(I,j,:), CS%uKoR(I,j,:), &
CS%uhEff(I,j,:), uFlx(I,j,:), &
CS%continuous_reconstruction, CS%remap_CS)
endif
enddo ; enddo
! y-flux
do J = G%jsc-1,G%jec ; do i = G%isc,G%iec
if (G%mask2dCv(i,J)>0.) then
call neutral_surface_flux(nk, CS%nsurf, CS%ppoly_deg, h(i,j,:), h(i,j+1,:), &
Tracer(i,j,:), Tracer(i,j+1,:), &
CS%vPoL(i,J,:), CS%vPoR(i,J,:), &
CS%vKoL(i,J,:), CS%vKoR(i,J,:), &
CS%vhEff(i,J,:), vFlx(i,J,:), &
CS%continuous_reconstruction, CS%remap_CS)
endif
enddo ; enddo
! Update the tracer concentration from divergence of neutral diffusive flux components
do j = G%jsc,G%jec ; do i = G%isc,G%iec
if (G%mask2dT(i,j)>0.) then
dTracer(:) = 0.
do ks = 1,CS%nsurf-1 ;
k = CS%uKoL(I,j,ks)
dTracer(k) = dTracer(k) + Coef_x(I,j) * uFlx(I,j,ks)
k = CS%uKoR(I-1,j,ks)
dTracer(k) = dTracer(k) - Coef_x(I-1,j) * uFlx(I-1,j,ks)
k = CS%vKoL(i,J,ks)
dTracer(k) = dTracer(k) + Coef_y(i,J) * vFlx(i,J,ks)
k = CS%vKoR(i,J-1,ks)
dTracer(k) = dTracer(k) - Coef_y(i,J-1) * vFlx(i,J-1,ks)
enddo
do k = 1, GV%ke
Tracer(i,j,k) = Tracer(i,j,k) + dTracer(k) * &
( G%IareaT(i,j) / ( h(i,j,k) + GV%H_subroundoff ) )
enddo
if(CS%id_neutral_diff_tracer_conc_tend(m) > 0 .or. &
CS%id_neutral_diff_tracer_cont_tend(m) > 0 .or. &
CS%id_neutral_diff_tracer_cont_tend_2d(m) > 0 ) then
do k = 1, GV%ke
tendency(i,j,k) = dTracer(k) * G%IareaT(i,j) * Idt
enddo
endif
endif
enddo ; enddo
! Diagnose vertically summed zonal flux, giving zonal tracer transport from ndiff.
! Note sign corresponds to downgradient flux convention.
if(CS%id_neutral_diff_tracer_trans_x_2d(m) > 0) then
do j = G%jsc,G%jec ; do I = G%isc-1,G%iec
trans_x_2d(I,j) = 0.
if (G%mask2dCu(I,j)>0.) then
do ks = 1,CS%nsurf-1 ;
trans_x_2d(I,j) = trans_x_2d(I,j) - Coef_x(I,j) * uFlx(I,j,ks)
enddo
trans_x_2d(I,j) = trans_x_2d(I,j) * Idt * convert
endif
enddo ; enddo
call post_data(CS%id_neutral_diff_tracer_trans_x_2d(m), trans_x_2d(:,:), CS%diag)
endif
! Diagnose vertically summed merid flux, giving meridional tracer transport from ndiff.
! Note sign corresponds to downgradient flux convention.
if(CS%id_neutral_diff_tracer_trans_y_2d(m) > 0) then
do J = G%jsc-1,G%jec ; do i = G%isc,G%iec
trans_y_2d(i,J) = 0.
if (G%mask2dCv(i,J)>0.) then
do ks = 1,CS%nsurf-1 ;
trans_y_2d(i,J) = trans_y_2d(i,J) - Coef_y(i,J) * vFlx(i,J,ks)
enddo
trans_y_2d(i,J) = trans_y_2d(i,J) * Idt * convert
endif
enddo ; enddo
call post_data(CS%id_neutral_diff_tracer_trans_y_2d(m), trans_y_2d(:,:), CS%diag)
endif
! post tendency of tracer content
if(CS%id_neutral_diff_tracer_cont_tend(m) > 0) then
call post_data(CS%id_neutral_diff_tracer_cont_tend(m), tendency(:,:,:)*convert, CS%diag)
endif
! post depth summed tendency for tracer content
if(CS%id_neutral_diff_tracer_cont_tend_2d(m) > 0) then
do j = G%jsc,G%jec ; do i = G%isc,G%iec
do k = 1, GV%ke
tendency_2d(i,j) = tendency_2d(i,j) + tendency(i,j,k)
enddo
enddo ; enddo
call post_data(CS%id_neutral_diff_tracer_cont_tend_2d(m), tendency_2d(:,:)*convert, CS%diag)
endif
! post tendency of tracer concentration; this step must be
! done after posting tracer content tendency, since we alter
! the tendency array.
if(CS%id_neutral_diff_tracer_conc_tend(m) > 0) then
do k = 1, GV%ke ; do j = G%jsc,G%jec ; do i = G%isc,G%iec
tendency(i,j,k) = tendency(i,j,k) / ( h(i,j,k) + GV%H_subroundoff )
enddo ; enddo ; enddo
call post_data(CS%id_neutral_diff_tracer_conc_tend(m), tendency, CS%diag)
endif
end subroutine neutral_diffusion
!> Returns interface scalar, Si, for a column of layer values, S.
subroutine interface_scalar(nk, h, S, Si, i_method)
integer, intent(in) :: nk !< Number of levels
real, dimension(nk), intent(in) :: h !< Layer thickness (H units)
real, dimension(nk), intent(in) :: S !< Layer scalar (conc, e.g. ppt)
real, dimension(nk+1), intent(inout) :: Si !< Interface scalar (conc, e.g. ppt)
integer, intent(in) :: i_method !< =1 use average of PLM edges
!! =2 use continuous PPM edge interpolation
! Local variables
integer :: k, km2, kp1
real, dimension(nk) :: diff
real :: Sb, Sa
call PLM_diff(nk, h, S, 2, 1, diff)
Si(1) = S(1) - 0.5 * diff(1)
if (i_method==1) then
do k = 2, nk
! Average of the two edge values (will be bounded and,
! when slopes are unlimited, notionally second-order accurate)
Sa = S(k-1) + 0.5 * diff(k-1) ! Lower edge value of a PLM reconstruction for layer above
Sb = S(k) - 0.5 * diff(k) ! Upper edge value of a PLM reconstruction for layer below
Si(k) = 0.5 * ( Sa + Sb )
enddo
elseif (i_method==2) then
do k = 2, nk
! PPM quasi-fourth order interpolation for edge values following
! equation 1.6 in Colella & Woodward, 1984: JCP 54, 174-201.
km2 = max(1, k-2)
kp1 = min(nk, k+1)
Si(k) = ppm_edge(h(km2), h(k-1), h(k), h(kp1), S(k-1), S(k), diff(k-1), diff(k))
enddo
endif
Si(nk+1) = S(nk) + 0.5 * diff(nk)
end subroutine interface_scalar
!> Returns the PPM quasi-fourth order edge value at k+1/2 following
!! equation 1.6 in Colella & Woodward, 1984: JCP 54, 174-201.
real function ppm_edge(hkm1, hk, hkp1, hkp2, Ak, Akp1, Pk, Pkp1)
real, intent(in) :: hkm1 !< Width of cell k-1
real, intent(in) :: hk !< Width of cell k
real, intent(in) :: hkp1 !< Width of cell k+1
real, intent(in) :: hkp2 !< Width of cell k+2
real, intent(in) :: Ak !< Average scalar value of cell k
real, intent(in) :: Akp1 !< Average scalar value of cell k+1
real, intent(in) :: Pk !< PLM slope for cell k
real, intent(in) :: Pkp1 !< PLM slope for cell k+1
! Local variables
real :: R_hk_hkp1, R_2hk_hkp1, R_hk_2hkp1, f1, f2, f3, f4
real, parameter :: h_neglect = 1.e-30
R_hk_hkp1 = hk + hkp1
if (R_hk_hkp1 <= 0.) then
ppm_edge = 0.5 * ( Ak + Akp1 )
return
endif
R_hk_hkp1 = 1. / R_hk_hkp1
if (hk<hkp1) then
ppm_edge = Ak + ( hk * R_hk_hkp1 ) * ( Akp1 - Ak )
else
ppm_edge = Akp1 + ( hkp1 * R_hk_hkp1 ) * ( Ak - Akp1 )
endif
R_2hk_hkp1 = 1. / ( ( 2. * hk + hkp1 ) + h_neglect )
R_hk_2hkp1 = 1. / ( ( hk + 2. * hkp1 ) + h_neglect )
f1 = 1./ ( ( hk + hkp1) + ( hkm1 + hkp2 ) )
f2 = 2. * ( hkp1 * hk ) * R_hk_hkp1 * &
( ( hkm1 + hk ) * R_2hk_hkp1 - ( hkp2 + hkp1 ) * R_hk_2hkp1 )
f3 = hk * ( hkm1 + hk ) * R_2hk_hkp1
f4 = hkp1 * ( hkp1 + hkp2 ) * R_hk_2hkp1
ppm_edge = ppm_edge + f1 * ( f2 * ( Akp1 - Ak ) - ( f3 * Pkp1 - f4 * Pk ) )
end function ppm_edge
!> Returns the average of a PPM reconstruction between two
!! fractional positions.
real function ppm_ave(xL, xR, aL, aR, aMean)
real, intent(in) :: xL !< Fraction position of left bound (0,1)
real, intent(in) :: xR !< Fraction position of right bound (0,1)
real, intent(in) :: aL !< Left edge scalar value, at x=0
real, intent(in) :: aR !< Right edge scalar value, at x=1
real, intent(in) :: aMean !< Average scalar value of cell
! Local variables
real :: dx, xave, a6, a6o3
dx = xR - xL
xave = 0.5 * ( xR + xL )
a6o3 = 2. * aMean - ( aL + aR ) ! a6 / 3.
a6 = 3. * a6o3
if (dx<0.) then
stop 'ppm_ave: dx<0 should not happend!'
elseif (dx>1.) then
stop 'ppm_ave: dx>1 should not happend!'
elseif (dx==0.) then
ppm_ave = aL + ( aR - aL ) * xR + a6 * xR * ( 1. - xR )
else
ppm_ave = ( aL + xave * ( ( aR - aL ) + a6 ) ) - a6o3 * ( xR**2 + xR * xL + xL**2 )
endif
end function ppm_ave
!> A true signum function that returns either -abs(a), when x<0; or abs(a) when x>0; or 0 when x=0.
real function signum(a,x)
real, intent(in) :: a !< The magnitude argument
real, intent(in) :: x !< The sign (or zero) argument
signum = sign(a,x)
if (x==0.) signum = 0.
end function signum
!> Returns PLM slopes for a column where the slopes are the difference in value across each cell.
!! The limiting follows equation 1.8 in Colella & Woodward, 1984: JCP 54, 174-201.
subroutine PLM_diff(nk, h, S, c_method, b_method, diff)
integer, intent(in) :: nk !< Number of levels
real, dimension(nk), intent(in) :: h !< Layer thickness (H units)
real, dimension(nk), intent(in) :: S !< Layer salinity (conc, e.g. ppt)
integer, intent(in) :: c_method !< Method to use for the centered difference
integer, intent(in) :: b_method !< =1, use PCM in first/last cell, =2 uses linear extrapolation
real, dimension(nk), intent(inout) :: diff !< Scalar difference across layer (conc, e.g. ppt)
!! determined by the following values for c_method:
!! 1. Second order finite difference (not recommended)
!! 2. Second order finite volume (used in original PPM)
!! 3. Finite-volume weighted least squares linear fit
!! \todo The use of c_method to choose a scheme is inefficient
!! and should eventually be moved up the call tree.
! Local variables
integer :: k
real :: hkm1, hk, hkp1, Skm1, Sk, Skp1, diff_l, diff_r, diff_c
do k = 2, nk-1
hkm1 = h(k-1)
hk = h(k)
hkp1 = h(k+1)
if ( ( hkp1 + hk ) * ( hkm1 + hk ) > 0.) then
Skm1 = S(k-1)
Sk = S(k)
Skp1 = S(k+1)
if (c_method==1) then
! Simple centered diff (from White)
if ( hk + 0.5 * (hkm1 + hkp1) /= 0. ) then
diff_c = ( Skp1 - Skm1 ) * ( hk / ( hk + 0.5 * (hkm1 + hkp1) ) )
else
diff_c = 0.
endif
elseif (c_method==2) then
! Second order accurate centered FV slope (from Colella and Woodward, JCP 1984)
diff_c = fv_diff(hkm1, hk, hkp1, Skm1, Sk, Skp1)
elseif (c_method==3) then
! Second order accurate finite-volume least squares slope
diff_c = hk * fvlsq_slope(hkm1, hk, hkp1, Skm1, Sk, Skp1)
endif
! Limit centered slope by twice the side differenced slopes
diff_l = 2. * ( Sk - Skm1 )
diff_r = 2. * ( Skp1 - Sk )
if ( signum(1., diff_l) * signum(1., diff_r) <= 0. ) then
diff(k) = 0. ! PCM for local extrema
else
diff(k) = sign( min( abs(diff_l), abs(diff_c), abs(diff_r) ), diff_c )
endif
else
diff(k) = 0. ! PCM next to vanished layers
endif
enddo
if (b_method==1) then ! PCM for top and bottom layer
diff(1) = 0.
diff(nk) = 0.
elseif (b_method==2) then ! Linear extrapolation for top and bottom interfaces
diff(1) = ( S(2) - S(1) ) * 2. * ( h(1) / ( h(1) + h(2) ) )
diff(nk) = S(nk) - S(nk-1) * 2. * ( h(nk) / ( h(nk-1) + h(nk) ) )
endif
end subroutine PLM_diff
!> Returns the cell-centered second-order finite volume (unlimited PLM) slope
!! using three consecutive cell widths and average values. Slope is returned
!! as a difference across the central cell (i.e. units of scalar S).
!! Discretization follows equation 1.7 in Colella & Woodward, 1984: JCP 54, 174-201.
real function fv_diff(hkm1, hk, hkp1, Skm1, Sk, Skp1)
real, intent(in) :: hkm1 !< Left cell width
real, intent(in) :: hk !< Center cell width
real, intent(in) :: hkp1 !< Right cell width
real, intent(in) :: Skm1 !< Left cell average value
real, intent(in) :: Sk !< Center cell average value
real, intent(in) :: Skp1 !< Right cell average value
! Local variables
real :: h_sum, hp, hm
h_sum = ( hkm1 + hkp1 ) + hk
if (h_sum /= 0.) h_sum = 1./ h_sum
hm = hkm1 + hk
if (hm /= 0.) hm = 1./ hm
hp = hkp1 + hk
if (hp /= 0.) hp = 1./ hp
fv_diff = ( hk * h_sum ) * &
( ( 2. * hkm1 + hk ) * hp * ( Skp1 - Sk ) &
+ ( 2. * hkp1 + hk ) * hm * ( Sk - Skm1 ) )
end function fv_diff
!> Returns the cell-centered second-order weighted least squares slope
!! using three consecutive cell widths and average values. Slope is returned
!! as a gradient (i.e. units of scalar S over width units).
real function fvlsq_slope(hkm1, hk, hkp1, Skm1, Sk, Skp1)
real, intent(in) :: hkm1 !< Left cell width
real, intent(in) :: hk !< Center cell width
real, intent(in) :: hkp1 !< Right cell width
real, intent(in) :: Skm1 !< Left cell average value
real, intent(in) :: Sk !< Center cell average value
real, intent(in) :: Skp1 !< Right cell average value
! Local variables
real :: xkm1, xkp1
real :: h_sum, hx_sum, hxsq_sum, hxy_sum, hy_sum, det
xkm1 = -0.5 * ( hk + hkm1 )
xkp1 = 0.5 * ( hk + hkp1 )
h_sum = ( hkm1 + hkp1 ) + hk
hx_sum = hkm1*xkm1 + hkp1*xkp1
hxsq_sum = hkm1*(xkm1**2) + hkp1*(xkp1**2)
hxy_sum = hkm1*xkm1*Skm1 + hkp1*xkp1*Skp1
hy_sum = ( hkm1*Skm1 + hkp1*Skp1 ) + hk*Sk
det = h_sum * hxsq_sum - hx_sum**2
if (det /= 0.) then
!a = ( hxsq_sum * hy_sum - hx_sum*hxy_sum ) / det ! a would be mean of straight line fit
fvlsq_slope = ( h_sum * hxy_sum - hx_sum*hy_sum ) / det ! Gradient of straight line fit
else
fvlsq_slope = 0. ! Adcroft's reciprocal rule
endif
end function fvlsq_slope
!> Returns positions within left/right columns of combined interfaces using continuous reconstructions of T/S
subroutine find_neutral_surface_positions_continuous(nk, Pl, Tl, Sl, dRdTl, dRdSl, Pr, Tr, Sr, dRdTr, dRdSr, PoL, &
PoR, KoL, KoR, hEff)
integer, intent(in) :: nk !< Number of levels
real, dimension(nk+1), intent(in) :: Pl !< Left-column interface pressure (Pa)
real, dimension(nk+1), intent(in) :: Tl !< Left-column interface potential temperature (degC)
real, dimension(nk+1), intent(in) :: Sl !< Left-column interface salinity (ppt)
real, dimension(nk+1), intent(in) :: dRdTl !< Left-column dRho/dT (kg/m3/degC)
real, dimension(nk+1), intent(in) :: dRdSl !< Left-column dRho/dS (kg/m3/ppt)
real, dimension(nk+1), intent(in) :: Pr !< Right-column interface pressure (Pa)
real, dimension(nk+1), intent(in) :: Tr !< Right-column interface potential temperature (degC)
real, dimension(nk+1), intent(in) :: Sr !< Right-column interface salinity (ppt)
real, dimension(nk+1), intent(in) :: dRdTr !< Left-column dRho/dT (kg/m3/degC)
real, dimension(nk+1), intent(in) :: dRdSr !< Left-column dRho/dS (kg/m3/ppt)
real, dimension(2*nk+2), intent(inout) :: PoL !< Fractional position of neutral surface within layer KoL of left column
real, dimension(2*nk+2), intent(inout) :: PoR !< Fractional position of neutral surface within layer KoR of right column
integer, dimension(2*nk+2), intent(inout) :: KoL !< Index of first left interface above neutral surface
integer, dimension(2*nk+2), intent(inout) :: KoR !< Index of first right interface above neutral surface
real, dimension(2*nk+1), intent(inout) :: hEff !< Effective thickness between two neutral surfaces (Pa)
! Local variables
integer :: ns ! Number of neutral surfaces
integer :: k_surface ! Index of neutral surface
integer :: kl ! Index of left interface
integer :: kr ! Index of right interface
real :: dRdT, dRdS ! dRho/dT and dRho/dS for the neutral surface
logical :: searching_left_column ! True if searching for the position of a right interface in the left column
logical :: searching_right_column ! True if searching for the position of a left interface in the right column
logical :: reached_bottom ! True if one of the bottom-most interfaces has been used as the target
integer :: krm1, klm1
real :: dRho, dRhoTop, dRhoBot, hL, hR
integer :: lastK_left, lastK_right
real :: lastP_left, lastP_right
ns = 2*nk+2
! Initialize variables for the search
kr = 1 ; lastK_right = 1 ; lastP_right = 0.
kl = 1 ; lastK_left = 1 ; lastP_left = 0.
reached_bottom = .false.
! Loop over each neutral surface, working from top to bottom
neutral_surfaces: do k_surface = 1, ns
klm1 = max(kl-1, 1)
if (klm1>nk) stop 'find_neutral_surface_positions(): klm1 went out of bounds!'
krm1 = max(kr-1, 1)
if (krm1>nk) stop 'find_neutral_surface_positions(): krm1 went out of bounds!'
! Potential density difference, rho(kr) - rho(kl)
dRho = 0.5 * ( ( dRdTr(kr) + dRdTl(kl) ) * ( Tr(kr) - Tl(kl) ) &
+ ( dRdSr(kr) + dRdSl(kl) ) * ( Sr(kr) - Sl(kl) ) )
! Which column has the lighter surface for the current indexes, kr and kl
if (.not. reached_bottom) then
if (dRho < 0.) then
searching_left_column = .true.
searching_right_column = .false.
elseif (dRho > 0.) then
searching_right_column = .true.
searching_left_column = .false.
else ! dRho == 0.
if (kl + kr == 2) then ! Still at surface
searching_left_column = .true.
searching_right_column = .false.
else ! Not the surface so we simply change direction
searching_left_column = .not. searching_left_column
searching_right_column = .not. searching_right_column
endif
endif
endif
if (searching_left_column) then
! Interpolate for the neutral surface position within the left column, layer klm1
! Potential density difference, rho(kl-1) - rho(kr) (should be negative)
dRhoTop = 0.5 * ( ( dRdTl(klm1) + dRdTr(kr) ) * ( Tl(klm1) - Tr(kr) ) &
+ ( dRdSl(klm1) + dRdSr(kr) ) * ( Sl(klm1) - Sr(kr) ) )
! Potential density difference, rho(kl) - rho(kr) (will be positive)
dRhoBot = 0.5 * ( ( dRdTl(klm1+1) + dRdTr(kr) ) * ( Tl(klm1+1) - Tr(kr) ) &
+ ( dRdSl(klm1+1) + dRdSr(kr) ) * ( Sl(klm1+1) - Sr(kr) ) )
! Because we are looking left, the right surface, kr, is lighter than klm1+1 and should be denser than klm1
! unless we are still at the top of the left column (kl=1)
if (dRhoTop > 0. .or. kr+kl==2) then
PoL(k_surface) = 0. ! The right surface is lighter than anything in layer klm1
elseif (dRhoTop >= dRhoBot) then ! Left layer is unstratified
PoL(k_surface) = 1.
else
! Linearly interpolate for the position between Pl(kl-1) and Pl(kl) where the density difference
! between right and left is zero.
PoL(k_surface) = interpolate_for_nondim_position( dRhoTop, Pl(klm1), dRhoBot, Pl(klm1+1) )
endif
if (PoL(k_surface)>=1. .and. klm1<nk) then ! >= is really ==, when PoL==1 we point to the bottom of the cell
klm1 = klm1 + 1
PoL(k_surface) = PoL(k_surface) - 1.
endif
if (real(klm1-lastK_left)+(PoL(k_surface)-lastP_left)<0.) then
PoL(k_surface) = lastP_left
klm1 = lastK_left
endif
KoL(k_surface) = klm1
if (kr <= nk) then
PoR(k_surface) = 0.
KoR(k_surface) = kr
else
PoR(k_surface) = 1.
KoR(k_surface) = nk
endif
if (kr <= nk) then
kr = kr + 1
else
reached_bottom = .true.
searching_right_column = .true.
searching_left_column = .false.
endif
elseif (searching_right_column) then
! Interpolate for the neutral surface position within the right column, layer krm1
! Potential density difference, rho(kr-1) - rho(kl) (should be negative)
dRhoTop = 0.5 * ( ( dRdTr(krm1) + dRdTl(kl) ) * ( Tr(krm1) - Tl(kl) ) &
+ ( dRdSr(krm1) + dRdSl(kl) ) * ( Sr(krm1) - Sl(kl) ) )
! Potential density difference, rho(kr) - rho(kl) (will be positive)
dRhoBot = 0.5 * ( ( dRdTr(krm1+1) + dRdTl(kl) ) * ( Tr(krm1+1) - Tl(kl) ) &
+ ( dRdSr(krm1+1) + dRdSl(kl) ) * ( Sr(krm1+1) - Sl(kl) ) )
! Because we are looking right, the left surface, kl, is lighter than krm1+1 and should be denser than krm1
! unless we are still at the top of the right column (kr=1)
if (dRhoTop >= 0. .or. kr+kl==2) then
PoR(k_surface) = 0. ! The left surface is lighter than anything in layer krm1
elseif (dRhoTop >= dRhoBot) then ! Right layer is unstratified
PoR(k_surface) = 1.
else
! Linearly interpolate for the position between Pr(kr-1) and Pr(kr) where the density difference
! between right and left is zero.
PoR(k_surface) = interpolate_for_nondim_position( dRhoTop, Pr(krm1), dRhoBot, Pr(krm1+1) )
endif