Newer
Older
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 icedyn_adv_umx
!!==============================================================================
!! *** MODULE icedyn_adv_umx ***
!! sea-ice : advection using the ULTIMATE-MACHO scheme
!!==============================================================================
!! History : 3.6 ! 2014-11 (C. Rousset, G. Madec) Original code
!! 4.0 ! 2018 (many people) SI3 [aka Sea Ice cube]
!!----------------------------------------------------------------------
#if defined key_si3
!!----------------------------------------------------------------------
!! 'key_si3' SI3 sea-ice model
!!----------------------------------------------------------------------
!! ice_dyn_adv_umx : update the tracer fields
!! ultimate_x(_y) : compute a tracer value at velocity points using ULTIMATE scheme at various orders
!! macho : compute the fluxes
!! nonosc_ice : limit the fluxes using a non-oscillatory algorithm
!!----------------------------------------------------------------------
USE phycst ! physical constant
USE dom_oce ! ocean domain
USE sbc_oce , ONLY : nn_fsbc ! update frequency of surface boundary condition
USE ice ! sea-ice variables
USE icevar ! sea-ice: operations
!
USE in_out_manager ! I/O manager
USE iom ! I/O manager library
USE lib_mpp ! MPP library
USE lib_fortran ! fortran utilities (glob_sum + no signed zero)
USE lbclnk ! lateral boundary conditions (or mpp links)
IMPLICIT NONE
PRIVATE
PUBLIC ice_dyn_adv_umx ! called by icedyn_adv.F90
!
INTEGER, PARAMETER :: np_advS = 1 ! advection for S and T: dVS/dt = -div( uVS ) => np_advS = 1
! or dVS/dt = -div( uA * uHS / u ) => np_advS = 2
! or dVS/dt = -div( uV * uS / u ) => np_advS = 3
INTEGER, PARAMETER :: np_limiter = 1 ! limiter: 1 = nonosc
! 2 = superbee
! 3 = h3
LOGICAL :: ll_upsxy = .TRUE. ! alternate directions for upstream
LOGICAL :: ll_hoxy = .TRUE. ! alternate directions for high order
LOGICAL :: ll_neg = .TRUE. ! if T interpolated at u/v points is negative or v_i < 1.e-6
! then interpolate T at u/v points using the upstream scheme
LOGICAL :: ll_prelim = .FALSE. ! prelimiter from: Zalesak(1979) eq. 14 => not well defined in 2D
!
REAL(wp) :: z1_6 = 1._wp / 6._wp ! =1/6
REAL(wp) :: z1_120 = 1._wp / 120._wp ! =1/120
!
INTEGER, ALLOCATABLE, DIMENSION(:,:,:) :: imsk_small, jmsk_small
!
!! * Substitutions
# include "do_loop_substitute.h90"
!!----------------------------------------------------------------------
!! NEMO/ICE 4.0 , NEMO Consortium (2018)
!! $Id: icedyn_adv_umx.F90 15049 2021-06-23 16:17:30Z clem $
!! Software governed by the CeCILL licence (./LICENSE)
!!----------------------------------------------------------------------
CONTAINS
SUBROUTINE ice_dyn_adv_umx( kn_umx, kt, pu_ice, pv_ice, ph_i, ph_s, ph_ip, &
& pato_i, pv_i, pv_s, psv_i, poa_i, pa_i, pa_ip, pv_ip, pv_il, pe_s, pe_i )
!!----------------------------------------------------------------------
!! *** ROUTINE ice_dyn_adv_umx ***
!!
!! ** Purpose : Compute the now trend due to total advection of
!! tracers and add it to the general trend of tracer equations
!! using an "Ultimate-Macho" scheme
!!
!! Reference : Leonard, B.P., 1991, Comput. Methods Appl. Mech. Eng., 88, 17-74.
!!----------------------------------------------------------------------
INTEGER , INTENT(in ) :: kn_umx ! order of the scheme (1-5=UM or 20=CEN2)
INTEGER , INTENT(in ) :: kt ! time step
REAL(wp), DIMENSION(:,:) , INTENT(in ) :: pu_ice ! ice i-velocity
REAL(wp), DIMENSION(:,:) , INTENT(in ) :: pv_ice ! ice j-velocity
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: ph_i ! ice thickness
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: ph_s ! snw thickness
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: ph_ip ! ice pond thickness
REAL(wp), DIMENSION(:,:) , INTENT(inout) :: pato_i ! open water area
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pv_i ! ice volume
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pv_s ! snw volume
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: psv_i ! salt content
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: poa_i ! age content
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pa_i ! ice concentration
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pa_ip ! melt pond concentration
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pv_ip ! melt pond volume
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pv_il ! melt pond lid volume
REAL(wp), DIMENSION(:,:,:,:), INTENT(inout) :: pe_s ! snw heat content
REAL(wp), DIMENSION(:,:,:,:), INTENT(inout) :: pe_i ! ice heat content
!
INTEGER :: ji, jj, jk, jl, jt ! dummy loop indices
INTEGER :: icycle ! number of sub-timestep for the advection
REAL(wp) :: zamsk ! 1 if advection of concentration, 0 if advection of other tracers
REAL(wp) :: zdt, z1_dt, zvi_cen
REAL(wp), DIMENSION(1) :: zcflprv, zcflnow ! for global communication
REAL(wp), DIMENSION(jpi,jpj) :: zudy, zvdx, zcu_box, zcv_box
REAL(wp), DIMENSION(jpi,jpj) :: zati1, zati2
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zu_cat, zv_cat
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zua_ho, zva_ho, zua_ups, zva_ups
REAL(wp), DIMENSION(jpi,jpj,jpl) :: z1_ai , z1_aip, zhvar
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zhi_max, zhs_max, zhip_max, zs_i, zsi_max
REAL(wp), DIMENSION(jpi,jpj,nlay_i,jpl) :: ze_i, zei_max
REAL(wp), DIMENSION(jpi,jpj,nlay_s,jpl) :: ze_s, zes_max
!
REAL(wp), ALLOCATABLE, DIMENSION(:,:,:) :: zuv_ho, zvv_ho, zuv_ups, zvv_ups, z1_vi, z1_vs
!! diagnostics
REAL(wp), DIMENSION(jpi,jpj) :: zdiag_adv_mass, zdiag_adv_salt, zdiag_adv_heat
!!----------------------------------------------------------------------
!
IF( kt == nit000 .AND. lwp ) WRITE(numout,*) '-- ice_dyn_adv_umx: Ultimate-Macho advection scheme'
!
! --- Record max of the surrounding 9-pts (for call Hbig) --- !
! thickness and salinity
WHERE( pv_i(:,:,:) >= epsi10 ) ; zs_i(:,:,:) = psv_i(:,:,:) / pv_i(:,:,:)
ELSEWHERE ; zs_i(:,:,:) = 0._wp
END WHERE
CALL icemax3D( ph_i , zhi_max )
CALL icemax3D( ph_s , zhs_max )
CALL icemax3D( ph_ip, zhip_max)
CALL icemax3D( zs_i , zsi_max )
CALL lbc_lnk( 'icedyn_adv_umx', zhi_max, 'T', 1._wp, zhs_max, 'T', 1._wp, zhip_max, 'T', 1._wp, zsi_max, 'T', 1._wp )
!
! enthalpies
DO jk = 1, nlay_i
WHERE( pv_i(:,:,:) >= epsi10 ) ; ze_i(:,:,jk,:) = pe_i(:,:,jk,:) / pv_i(:,:,:)
ELSEWHERE ; ze_i(:,:,jk,:) = 0._wp
END WHERE
END DO
DO jk = 1, nlay_s
WHERE( pv_s(:,:,:) >= epsi10 ) ; ze_s(:,:,jk,:) = pe_s(:,:,jk,:) / pv_s(:,:,:)
ELSEWHERE ; ze_s(:,:,jk,:) = 0._wp
END WHERE
END DO
CALL icemax4D( ze_i , zei_max )
CALL icemax4D( ze_s , zes_max )
CALL lbc_lnk( 'icedyn_adv_umx', zei_max, 'T', 1._wp )
CALL lbc_lnk( 'icedyn_adv_umx', zes_max, 'T', 1._wp )
!
!
! --- If ice drift is too fast, use subtime steps for advection (CFL test for stability) --- !
! Note: the advection split is applied at the next time-step in order to avoid blocking global comm.
! this should not affect too much the stability
zcflnow(1) = MAXVAL( ABS( pu_ice(:,:) ) * rDt_ice * r1_e1u(:,:) )
zcflnow(1) = MAX( zcflnow(1), MAXVAL( ABS( pv_ice(:,:) ) * rDt_ice * r1_e2v(:,:) ) )
! non-blocking global communication send zcflnow and receive zcflprv
CALL mpp_delay_max( 'icedyn_adv_umx', 'cflice', zcflnow(:), zcflprv(:), kt == nitend - nn_fsbc + 1 )
IF( zcflprv(1) > .5 ) THEN ; icycle = 2
ELSE ; icycle = 1
ENDIF
zdt = rDt_ice / REAL(icycle)
z1_dt = 1._wp / zdt
! --- transport --- !
zudy(:,:) = pu_ice(:,:) * e2u(:,:)
zvdx(:,:) = pv_ice(:,:) * e1v(:,:)
!
! setup transport for each ice cat
DO jl = 1, jpl
zu_cat(:,:,jl) = zudy(:,:)
zv_cat(:,:,jl) = zvdx(:,:)
END DO
!
! --- define velocity for advection: u*grad(H) --- !
DO_2D( nn_hls-1, nn_hls, nn_hls, nn_hls )
IF ( pu_ice(ji,jj) * pu_ice(ji-1,jj) <= 0._wp ) THEN ; zcu_box(ji,jj) = 0._wp
ELSEIF( pu_ice(ji,jj) > 0._wp ) THEN ; zcu_box(ji,jj) = pu_ice(ji-1,jj)
ELSE ; zcu_box(ji,jj) = pu_ice(ji ,jj)
ENDIF
END_2D
DO_2D( nn_hls, nn_hls, nn_hls-1, nn_hls )
IF ( pv_ice(ji,jj) * pv_ice(ji,jj-1) <= 0._wp ) THEN ; zcv_box(ji,jj) = 0._wp
ELSEIF( pv_ice(ji,jj) > 0._wp ) THEN ; zcv_box(ji,jj) = pv_ice(ji,jj-1)
ELSE ; zcv_box(ji,jj) = pv_ice(ji,jj )
ENDIF
END_2D
!---------------!
!== advection ==!
!---------------!
DO jt = 1, icycle
! diagnostics
zdiag_adv_mass(:,:) = SUM( pv_i (:,:,:) , dim=3 ) * rhoi + SUM( pv_s (:,:,:) , dim=3 ) * rhos &
& + SUM( pv_ip(:,:,:) , dim=3 ) * rhow + SUM( pv_il(:,:,:) , dim=3 ) * rhow
zdiag_adv_salt(:,:) = SUM( psv_i(:,:,:) , dim=3 ) * rhoi
zdiag_adv_heat(:,:) = - SUM(SUM( pe_i(:,:,1:nlay_i,:) , dim=4 ), dim=3 ) &
& - SUM(SUM( pe_s(:,:,1:nlay_s,:) , dim=4 ), dim=3 )
! record at_i before advection (for open water)
zati1(:,:) = SUM( pa_i(:,:,:), dim=3 )
! inverse of A and Ap
WHERE( pa_i(:,:,:) >= epsi20 ) ; z1_ai(:,:,:) = 1._wp / pa_i(:,:,:)
ELSEWHERE ; z1_ai(:,:,:) = 0.
END WHERE
WHERE( pa_ip(:,:,:) >= epsi20 ) ; z1_aip(:,:,:) = 1._wp / pa_ip(:,:,:)
ELSEWHERE ; z1_aip(:,:,:) = 0.
END WHERE
!
! setup a mask where advection will be upstream
IF( ll_neg ) THEN
IF( .NOT. ALLOCATED(imsk_small) ) ALLOCATE( imsk_small(jpi,jpj,jpl) )
IF( .NOT. ALLOCATED(jmsk_small) ) ALLOCATE( jmsk_small(jpi,jpj,jpl) )
DO jl = 1, jpl
DO_2D( 1, 0, nn_hls, nn_hls )
zvi_cen = 0.5_wp * ( pv_i(ji+1,jj,jl) + pv_i(ji,jj,jl) )
IF( zvi_cen < epsi06) THEN ; imsk_small(ji,jj,jl) = 0
ELSE ; imsk_small(ji,jj,jl) = 1 ; ENDIF
END_2D
DO_2D( nn_hls, nn_hls, 1, 0 )
zvi_cen = 0.5_wp * ( pv_i(ji,jj+1,jl) + pv_i(ji,jj,jl) )
IF( zvi_cen < epsi06) THEN ; jmsk_small(ji,jj,jl) = 0
ELSE ; jmsk_small(ji,jj,jl) = 1 ; ENDIF
END_2D
END DO
ENDIF
!
! ----------------------- !
! ==> start advection <== !
! ----------------------- !
!
!== Ice area ==!
zamsk = 1._wp
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy, zvdx, zu_cat , zv_cat , zcu_box, zcv_box, &
& pa_i, pa_i, zua_ups, zva_ups, zua_ho , zva_ho )
!
! ! --------------------------------- !
IF( np_advS == 1 ) THEN ! -- advection form: -div( uVS ) -- !
! ! --------------------------------- !
zamsk = 0._wp
!== Ice volume ==!
zhvar(:,:,:) = pv_i(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_i, zua_ups, zva_ups )
!== Snw volume ==!
zhvar(:,:,:) = pv_s(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_s, zua_ups, zva_ups )
!
zamsk = 1._wp
!== Salt content ==!
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zu_cat, zv_cat, zcu_box, zcv_box, &
& psv_i, psv_i )
!== Ice heat content ==!
DO jk = 1, nlay_i
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zu_cat, zv_cat, zcu_box, zcv_box, &
& pe_i(:,:,jk,:), pe_i(:,:,jk,:) )
END DO
!== Snw heat content ==!
DO jk = 1, nlay_s
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zu_cat, zv_cat, zcu_box, zcv_box, &
& pe_s(:,:,jk,:), pe_s(:,:,jk,:) )
END DO
!
! ! ------------------------------------------ !
ELSEIF( np_advS == 2 ) THEN ! -- advection form: -div( uA * uHS / u ) -- !
! ! ------------------------------------------ !
zamsk = 0._wp
!== Ice volume ==!
zhvar(:,:,:) = pv_i(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_i, zua_ups, zva_ups )
!== Snw volume ==!
zhvar(:,:,:) = pv_s(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_s, zua_ups, zva_ups )
!== Salt content ==!
zhvar(:,:,:) = psv_i(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, psv_i, zua_ups, zva_ups )
!== Ice heat content ==!
DO jk = 1, nlay_i
zhvar(:,:,:) = pe_i(:,:,jk,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho, zva_ho, zcu_box, zcv_box, &
& zhvar, pe_i(:,:,jk,:), zua_ups, zva_ups )
END DO
!== Snw heat content ==!
DO jk = 1, nlay_s
zhvar(:,:,:) = pe_s(:,:,jk,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho, zva_ho, zcu_box, zcv_box, &
& zhvar, pe_s(:,:,jk,:), zua_ups, zva_ups )
END DO
!
! ! ----------------------------------------- !
ELSEIF( np_advS == 3 ) THEN ! -- advection form: -div( uV * uS / u ) -- !
! ! ----------------------------------------- !
zamsk = 0._wp
!
ALLOCATE( zuv_ho (jpi,jpj,jpl), zvv_ho (jpi,jpj,jpl), &
& zuv_ups(jpi,jpj,jpl), zvv_ups(jpi,jpj,jpl), z1_vi(jpi,jpj,jpl), z1_vs(jpi,jpj,jpl) )
!
! inverse of Vi
WHERE( pv_i(:,:,:) >= epsi20 ) ; z1_vi(:,:,:) = 1._wp / pv_i(:,:,:)
ELSEWHERE ; z1_vi(:,:,:) = 0.
END WHERE
! inverse of Vs
WHERE( pv_s(:,:,:) >= epsi20 ) ; z1_vs(:,:,:) = 1._wp / pv_s(:,:,:)
ELSEWHERE ; z1_vs(:,:,:) = 0.
END WHERE
!
! It is important to first calculate the ice fields and then the snow fields (because we use the same arrays)
!
!== Ice volume ==!
zuv_ups = zua_ups
zvv_ups = zva_ups
zhvar(:,:,:) = pv_i(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_i, zuv_ups, zvv_ups, zuv_ho , zvv_ho )
!== Salt content ==!
zhvar(:,:,:) = psv_i(:,:,:) * z1_vi(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zuv_ho , zvv_ho , zcu_box, zcv_box, &
& zhvar, psv_i, zuv_ups, zvv_ups )
!== Ice heat content ==!
DO jk = 1, nlay_i
zhvar(:,:,:) = pe_i(:,:,jk,:) * z1_vi(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zuv_ho, zvv_ho, zcu_box, zcv_box, &
& zhvar, pe_i(:,:,jk,:), zuv_ups, zvv_ups )
END DO
!== Snow volume ==!
zuv_ups = zua_ups
zvv_ups = zva_ups
zhvar(:,:,:) = pv_s(:,:,:) * z1_ai(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_s, zuv_ups, zvv_ups, zuv_ho , zvv_ho )
!== Snw heat content ==!
DO jk = 1, nlay_s
zhvar(:,:,:) = pe_s(:,:,jk,:) * z1_vs(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx, zuv_ho, zvv_ho, zcu_box, zcv_box, &
& zhvar, pe_s(:,:,jk,:), zuv_ups, zvv_ups )
END DO
!
DEALLOCATE( zuv_ho, zvv_ho, zuv_ups, zvv_ups, z1_vi, z1_vs )
!
ENDIF
!
!== Ice age ==!
zamsk = 1._wp
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zu_cat, zv_cat, zcu_box, zcv_box, &
& poa_i, poa_i )
!
!== melt ponds ==!
IF ( ln_pnd_LEV .OR. ln_pnd_TOPO ) THEN
! concentration
zamsk = 1._wp
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zu_cat , zv_cat , zcu_box, zcv_box, &
& pa_ip, pa_ip, zua_ups, zva_ups, zua_ho , zva_ho )
! volume
zamsk = 0._wp
zhvar(:,:,:) = pv_ip(:,:,:) * z1_aip(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_ip, zua_ups, zva_ups )
! lid
IF ( ln_pnd_lids ) THEN
zamsk = 0._wp
zhvar(:,:,:) = pv_il(:,:,:) * z1_aip(:,:,:)
CALL adv_umx( zamsk, kn_umx, jt, kt, zdt, zudy , zvdx , zua_ho , zva_ho , zcu_box, zcv_box, &
& zhvar, pv_il, zua_ups, zva_ups )
ENDIF
ENDIF
! --- Lateral boundary conditions --- !
IF ( ( ln_pnd_LEV .OR. ln_pnd_TOPO ) .AND. ln_pnd_lids ) THEN
CALL lbc_lnk( 'icedyn_adv_umx', pa_i,'T',1._wp, pv_i,'T',1._wp, pv_s,'T',1._wp, psv_i,'T',1._wp, poa_i,'T',1._wp &
& , pa_ip,'T',1._wp, pv_ip,'T',1._wp, pv_il,'T',1._wp )
ELSEIF( ( ln_pnd_LEV .OR. ln_pnd_TOPO ) .AND. .NOT.ln_pnd_lids ) THEN
CALL lbc_lnk( 'icedyn_adv_umx', pa_i,'T',1._wp, pv_i,'T',1._wp, pv_s,'T',1._wp, psv_i,'T',1._wp, poa_i,'T',1._wp &
& , pa_ip,'T',1._wp, pv_ip,'T',1._wp )
ELSE
CALL lbc_lnk( 'icedyn_adv_umx', pa_i,'T',1._wp, pv_i,'T',1._wp, pv_s,'T',1._wp, psv_i,'T',1._wp, poa_i,'T',1._wp )
ENDIF
CALL lbc_lnk( 'icedyn_adv_umx', pe_i, 'T', 1._wp )
CALL lbc_lnk( 'icedyn_adv_umx', pe_s, 'T', 1._wp )
!
!== Open water area ==!
zati2(:,:) = SUM( pa_i(:,:,:), dim=3 )
DO_2D( 0, 0, 0, 0 )
pato_i(ji,jj) = pato_i(ji,jj) - ( zati2(ji,jj) - zati1(ji,jj) ) &
& - ( zudy(ji,jj) - zudy(ji-1,jj) + zvdx(ji,jj) - zvdx(ji,jj-1) ) * r1_e1e2t(ji,jj) * zdt
END_2D
CALL lbc_lnk( 'icedyn_adv_umx', pato_i, 'T', 1._wp )
!
! --- diagnostics --- !
diag_adv_mass(:,:) = diag_adv_mass(:,:) + ( SUM( pv_i (:,:,:) , dim=3 ) * rhoi + SUM( pv_s (:,:,:) , dim=3 ) * rhos &
& + SUM( pv_ip(:,:,:) , dim=3 ) * rhow + SUM( pv_il(:,:,:) , dim=3 ) * rhow &
& - zdiag_adv_mass(:,:) ) * z1_dt
diag_adv_salt(:,:) = diag_adv_salt(:,:) + ( SUM( psv_i(:,:,:) , dim=3 ) * rhoi &
& - zdiag_adv_salt(:,:) ) * z1_dt
diag_adv_heat(:,:) = diag_adv_heat(:,:) + ( - SUM(SUM( pe_i(:,:,1:nlay_i,:) , dim=4 ), dim=3 ) &
& - SUM(SUM( pe_s(:,:,1:nlay_s,:) , dim=4 ), dim=3 ) &
& - zdiag_adv_heat(:,:) ) * z1_dt
!
! --- Ensure non-negative fields and in-bound thicknesses --- !
! Remove negative values (conservation is ensured)
! (because advected fields are not perfectly bounded and tiny negative values can occur, e.g. -1.e-20)
CALL ice_var_zapneg( zdt, pato_i, pv_i, pv_s, psv_i, poa_i, pa_i, pa_ip, pv_ip, pv_il, pe_s, pe_i )
!
! --- Make sure ice thickness is not too big --- !
! (because ice thickness can be too large where ice concentration is very small)
CALL Hbig( zdt, zhi_max, zhs_max, zhip_max, zsi_max, zes_max, zei_max, &
& pv_i, pv_s, pa_i, pa_ip, pv_ip, psv_i, pe_s, pe_i )
!
! --- Ensure snow load is not too big --- !
CALL Hsnow( zdt, pv_i, pv_s, pa_i, pa_ip, pe_s )
!
END DO
!
END SUBROUTINE ice_dyn_adv_umx
SUBROUTINE adv_umx( pamsk, kn_umx, jt, kt, pdt, pu, pv, puc, pvc, pubox, pvbox, &
& pt, ptc, pua_ups, pva_ups, pua_ho, pva_ho )
!!----------------------------------------------------------------------
!! *** ROUTINE adv_umx ***
!!
!! ** Purpose : Compute the now trend due to total advection of
!! tracers and add it to the general trend of tracer equations
!!
!! ** Method : - calculate upstream fluxes and upstream solution for tracers V/A(=H) etc
!! - calculate tracer H at u and v points (Ultimate)
!! - calculate the high order fluxes using alterning directions (Macho)
!! - apply a limiter on the fluxes (nonosc_ice)
!! - convert this tracer flux to a "volume" flux (uH -> uV)
!! - apply a limiter a second time on the volumes fluxes (nonosc_ice)
!! - calculate the high order solution for V
!!
!! ** Action : solve 3 equations => a) dA/dt = -div(uA)
!! b) dV/dt = -div(uV) using dH/dt = -u.grad(H)
!! c) dVS/dt = -div(uVS) using either dHS/dt = -u.grad(HS) or dS/dt = -u.grad(S)
!!
!! in eq. b), - fluxes uH are evaluated (with UMx) and limited with nonosc_ice. This step is necessary to get a good H.
!! - then we convert this flux to a "volume" flux this way => uH * uA / u
!! where uA is the flux from eq. a)
!! this "volume" flux is also limited with nonosc_ice (otherwise overshoots can occur)
!! - at last we estimate dV/dt = -div(uH * uA / u)
!!
!! in eq. c), one can solve the equation for S (ln_advS=T), then dVS/dt = -div(uV * uS / u)
!! or for HS (ln_advS=F), then dVS/dt = -div(uA * uHS / u)
!!
!! ** Note : - this method can lead to tiny negative V (-1.e-20) => set it to 0 while conserving mass etc.
!! - At the ice edge, Ultimate scheme can lead to:
!! 1) negative interpolated tracers at u-v points
!! 2) non-zero interpolated tracers at u-v points eventhough there is no ice and velocity is outward
!! Solution for 1): apply an upstream scheme when it occurs. A better solution would be to degrade the order of
!! the scheme automatically by applying a mask of the ice cover inside Ultimate (not done).
!! Solution for 2): we set it to 0 in this case
!! - Eventhough 1D tests give very good results (typically the one from Schar & Smolarkiewiecz), the 2D is less good.
!! Large values of H can appear for very small ice concentration, and when it does it messes the things up since we
!! work on H (and not V). It is partly related to the multi-category approach
!! Therefore, after advection we limit the thickness to the largest value of the 9-points around (only if ice
!! concentration is small). We also limit S and T.
!!----------------------------------------------------------------------
REAL(wp) , INTENT(in ) :: pamsk ! advection of concentration (1) or other tracers (0)
INTEGER , INTENT(in ) :: kn_umx ! order of the scheme (1-5=UM or 20=CEN2)
INTEGER , INTENT(in ) :: jt ! number of sub-iteration
INTEGER , INTENT(in ) :: kt ! number of iteration
REAL(wp) , INTENT(in ) :: pdt ! tracer time-step
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pu , pv ! 2 ice velocity components => u*e2
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: puc , pvc ! 2 ice velocity components => u*e2 or u*a*e2u
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pubox, pvbox ! upstream velocity
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: pt ! tracer field
REAL(wp), DIMENSION(:,:,:) , INTENT(inout) :: ptc ! tracer content field
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT(inout), OPTIONAL :: pua_ups, pva_ups ! upstream u*a fluxes
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out), OPTIONAL :: pua_ho, pva_ho ! high order u*a fluxes
!
INTEGER :: ji, jj, jl ! dummy loop indices
REAL(wp) :: ztra ! local scalar
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zfu_ho , zfv_ho , zpt
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zfu_ups, zfv_ups, zt_ups
!!----------------------------------------------------------------------
!
! Upstream (_ups) fluxes
! -----------------------
CALL upstream( pamsk, jt, kt, pdt, pt, pu, pv, zt_ups, zfu_ups, zfv_ups )
! High order (_ho) fluxes
! -----------------------
SELECT CASE( kn_umx )
!
CASE ( 20 ) !== centered second order ==!
!
CALL cen2( pamsk, jt, kt, pdt, pt, pu, pv, zt_ups, zfu_ups, zfv_ups, zfu_ho, zfv_ho )
!
CASE ( 1:5 ) !== 1st to 5th order ULTIMATE-MACHO scheme ==!
!
CALL macho( pamsk, kn_umx, jt, kt, pdt, pt, pu, pv, pubox, pvbox, zt_ups, zfu_ups, zfv_ups, zfu_ho, zfv_ho )
!
END SELECT
!
! --ho --ho
! new fluxes = u*H * u*a / u
! ----------------------------
IF( pamsk == 0._wp ) THEN
DO jl = 1, jpl
DO_2D( 1, 0, 0, 0 )
IF( ABS( pu(ji,jj) ) > epsi10 ) THEN
zfu_ho (ji,jj,jl) = zfu_ho (ji,jj,jl) * puc (ji,jj,jl) / pu(ji,jj)
zfu_ups(ji,jj,jl) = zfu_ups(ji,jj,jl) * pua_ups(ji,jj,jl) / pu(ji,jj)
ELSE
zfu_ho (ji,jj,jl) = 0._wp
zfu_ups(ji,jj,jl) = 0._wp
ENDIF
!
END_2D
DO_2D( 0, 0, 1, 0 )
IF( ABS( pv(ji,jj) ) > epsi10 ) THEN
zfv_ho (ji,jj,jl) = zfv_ho (ji,jj,jl) * pvc (ji,jj,jl) / pv(ji,jj)
zfv_ups(ji,jj,jl) = zfv_ups(ji,jj,jl) * pva_ups(ji,jj,jl) / pv(ji,jj)
ELSE
zfv_ho (ji,jj,jl) = 0._wp
zfv_ups(ji,jj,jl) = 0._wp
ENDIF
END_2D
END DO
! the new "volume" fluxes must also be "flux corrected"
! thus we calculate the upstream solution and apply a limiter again
DO jl = 1, jpl
DO_2D( 0, 0, 0, 0 )
ztra = - ( zfu_ups(ji,jj,jl) - zfu_ups(ji-1,jj,jl) + zfv_ups(ji,jj,jl) - zfv_ups(ji,jj-1,jl) )
!
zt_ups(ji,jj,jl) = ( ptc(ji,jj,jl) + ztra * r1_e1e2t(ji,jj) * pdt ) * tmask(ji,jj,1)
END_2D
END DO
CALL lbc_lnk( 'icedyn_adv_umx', zt_ups, 'T', 1.0_wp )
!
IF ( np_limiter == 1 ) THEN
CALL nonosc_ice( 1._wp, pdt, pu, pv, ptc, zt_ups, zfu_ups, zfv_ups, zfu_ho, zfv_ho )
ELSEIF( np_limiter == 2 .OR. np_limiter == 3 ) THEN
CALL limiter_x( pdt, pu, ptc, zfu_ups, zfu_ho )
CALL limiter_y( pdt, pv, ptc, zfv_ups, zfv_ho )
ENDIF
!
ENDIF
! --ho --ups
! in case of advection of A: output u*a and u*a
! -----------------------------------------------
IF( PRESENT( pua_ho ) ) THEN
DO jl = 1, jpl
DO_2D( 1, 0, 0, 0 )
pua_ho (ji,jj,jl) = zfu_ho (ji,jj,jl)
pua_ups(ji,jj,jl) = zfu_ups(ji,jj,jl)
END_2D
DO_2D( 0, 0, 1, 0 )
pva_ho (ji,jj,jl) = zfv_ho (ji,jj,jl)
pva_ups(ji,jj,jl) = zfv_ups(ji,jj,jl)
END_2D
END DO
ENDIF
!
! final trend with corrected fluxes
! ---------------------------------
DO jl = 1, jpl
DO_2D( 0, 0, 0, 0 )
ztra = - ( zfu_ho(ji,jj,jl) - zfu_ho(ji-1,jj,jl) + zfv_ho(ji,jj,jl) - zfv_ho(ji,jj-1,jl) )
!
ptc(ji,jj,jl) = ( ptc(ji,jj,jl) + ztra * r1_e1e2t(ji,jj) * pdt ) * tmask(ji,jj,1)
END_2D
END DO
!
END SUBROUTINE adv_umx
SUBROUTINE upstream( pamsk, jt, kt, pdt, pt, pu, pv, pt_ups, pfu_ups, pfv_ups )
!!---------------------------------------------------------------------
!! *** ROUTINE upstream ***
!!
!! ** Purpose : compute the upstream fluxes and upstream guess of tracer
!!----------------------------------------------------------------------
REAL(wp) , INTENT(in ) :: pamsk ! advection of concentration (1) or other tracers (0)
INTEGER , INTENT(in ) :: jt ! number of sub-iteration
INTEGER , INTENT(in ) :: kt ! number of iteration
REAL(wp) , INTENT(in ) :: pdt ! tracer time-step
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pt ! tracer fields
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pu, pv ! 2 ice velocity components
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out) :: pt_ups ! upstream guess of tracer
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out) :: pfu_ups, pfv_ups ! upstream fluxes
!
INTEGER :: ji, jj, jl ! dummy loop indices
REAL(wp) :: ztra ! local scalar
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zpt
!!----------------------------------------------------------------------
IF( .NOT. ll_upsxy ) THEN !** no alternate directions **!
!
DO jl = 1, jpl
DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls-1 )
pfu_ups(ji,jj,jl) = MAX( pu(ji,jj), 0._wp ) * pt(ji,jj,jl) + MIN( pu(ji,jj), 0._wp ) * pt(ji+1,jj,jl)
pfv_ups(ji,jj,jl) = MAX( pv(ji,jj), 0._wp ) * pt(ji,jj,jl) + MIN( pv(ji,jj), 0._wp ) * pt(ji,jj+1,jl)
END_2D
END DO
!
ELSE !** alternate directions **!
!
IF( MOD( (kt - 1) / nn_fsbc , 2 ) == MOD( (jt - 1) , 2 ) ) THEN !== odd ice time step: adv_x then adv_y ==!
!
DO jl = 1, jpl !-- flux in x-direction
DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls )
pfu_ups(ji,jj,jl) = MAX( pu(ji,jj), 0._wp ) * pt(ji,jj,jl) + MIN( pu(ji,jj), 0._wp ) * pt(ji+1,jj,jl)
END_2D
END DO
!
DO jl = 1, jpl !-- first guess of tracer from u-flux
DO_2D( nn_hls-1, nn_hls-1, nn_hls, nn_hls )
ztra = - ( pfu_ups(ji,jj,jl) - pfu_ups(ji-1,jj,jl) ) &
& + ( pu (ji,jj ) - pu (ji-1,jj ) ) * pt(ji,jj,jl) * (1.-pamsk)
!
zpt(ji,jj,jl) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1)
END_2D
END DO
!
DO jl = 1, jpl !-- flux in y-direction
DO_2D( nn_hls-1, nn_hls-1, nn_hls, nn_hls-1 )
pfv_ups(ji,jj,jl) = MAX( pv(ji,jj), 0._wp ) * zpt(ji,jj,jl) + MIN( pv(ji,jj), 0._wp ) * zpt(ji,jj+1,jl)
END_2D
END DO
!
ELSE !== even ice time step: adv_y then adv_x ==!
!
DO jl = 1, jpl !-- flux in y-direction
DO_2D( nn_hls, nn_hls, nn_hls, nn_hls-1 )
pfv_ups(ji,jj,jl) = MAX( pv(ji,jj), 0._wp ) * pt(ji,jj,jl) + MIN( pv(ji,jj), 0._wp ) * pt(ji,jj+1,jl)
END_2D
END DO
!
DO jl = 1, jpl !-- first guess of tracer from v-flux
DO_2D( nn_hls, nn_hls, nn_hls-1, nn_hls-1 )
ztra = - ( pfv_ups(ji,jj,jl) - pfv_ups(ji,jj-1,jl) ) &
& + ( pv (ji,jj ) - pv (ji,jj-1 ) ) * pt(ji,jj,jl) * (1.-pamsk)
!
zpt(ji,jj,jl) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1)
END_2D
END DO
!
DO jl = 1, jpl !-- flux in x-direction
DO_2D( nn_hls, nn_hls-1, nn_hls-1, nn_hls-1 )
pfu_ups(ji,jj,jl) = MAX( pu(ji,jj), 0._wp ) * zpt(ji,jj,jl) + MIN( pu(ji,jj), 0._wp ) * zpt(ji+1,jj,jl)
END_2D
END DO
!
ENDIF
ENDIF
!
DO jl = 1, jpl !-- after tracer with upstream scheme
DO_2D( nn_hls-1, nn_hls-1, nn_hls-1, nn_hls-1 )
ztra = - ( pfu_ups(ji,jj,jl) - pfu_ups(ji-1,jj ,jl) &
& + pfv_ups(ji,jj,jl) - pfv_ups(ji ,jj-1,jl) ) &
& + ( pu (ji,jj ) - pu (ji-1,jj ) &
& + pv (ji,jj ) - pv (ji ,jj-1 ) ) * pt(ji,jj,jl) * (1.-pamsk)
!
pt_ups(ji,jj,jl) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1)
END_2D
END DO
IF( nn_hls == 1 ) CALL lbc_lnk( 'icedyn_adv_umx', pt_ups, 'T', 1.0_wp )
END SUBROUTINE upstream
SUBROUTINE cen2( pamsk, jt, kt, pdt, pt, pu, pv, pt_ups, pfu_ups, pfv_ups, pfu_ho, pfv_ho )
!!---------------------------------------------------------------------
!! *** ROUTINE cen2 ***
!!
!! ** Purpose : compute the high order fluxes using a centered
!! second order scheme
!!----------------------------------------------------------------------
REAL(wp) , INTENT(in ) :: pamsk ! advection of concentration (1) or other tracers (0)
INTEGER , INTENT(in ) :: jt ! number of sub-iteration
INTEGER , INTENT(in ) :: kt ! number of iteration
REAL(wp) , INTENT(in ) :: pdt ! tracer time-step
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pt ! tracer fields
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pu, pv ! 2 ice velocity components
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pt_ups ! upstream guess of tracer
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pfu_ups, pfv_ups ! upstream fluxes
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out) :: pfu_ho, pfv_ho ! high order fluxes
!
INTEGER :: ji, jj, jl ! dummy loop indices
REAL(wp) :: ztra ! local scalar
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zpt
!!----------------------------------------------------------------------
!
IF( .NOT.ll_hoxy ) THEN !** no alternate directions **!
!
DO jl = 1, jpl
DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls )
pfu_ho(ji,jj,jl) = 0.5_wp * pu(ji,jj) * ( pt(ji,jj,jl) + pt(ji+1,jj ,jl) )
END_2D
DO_2D( nn_hls, nn_hls, nn_hls, nn_hls-1 )
pfv_ho(ji,jj,jl) = 0.5_wp * pv(ji,jj) * ( pt(ji,jj,jl) + pt(ji ,jj+1,jl) )
END_2D
END DO
!
IF ( np_limiter == 1 ) THEN
CALL nonosc_ice( pamsk, pdt, pu, pv, pt, pt_ups, pfu_ups, pfv_ups, pfu_ho, pfv_ho )
ELSEIF( np_limiter == 2 .OR. np_limiter == 3 ) THEN
CALL limiter_x( pdt, pu, pt, pfu_ups, pfu_ho )
CALL limiter_y( pdt, pv, pt, pfv_ups, pfv_ho )
ENDIF
!
ELSE !** alternate directions **!
!
IF( MOD( (kt - 1) / nn_fsbc , 2 ) == MOD( (jt - 1) , 2 ) ) THEN !== odd ice time step: adv_x then adv_y ==!
!
DO jl = 1, jpl !-- flux in x-direction
DO_2D( nn_hls, nn_hls-1, nn_hls, nn_hls )
pfu_ho(ji,jj,jl) = 0.5_wp * pu(ji,jj) * ( pt(ji,jj,jl) + pt(ji+1,jj,jl) )
END_2D
END DO
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_x( pdt, pu, pt, pfu_ups, pfu_ho )
DO jl = 1, jpl !-- first guess of tracer from u-flux
DO_2D( nn_hls-1, nn_hls-1, nn_hls, nn_hls )
ztra = - ( pfu_ho(ji,jj,jl) - pfu_ho(ji-1,jj,jl) ) &
& + ( pu (ji,jj ) - pu (ji-1,jj ) ) * pt(ji,jj,jl) * (1.-pamsk)
!
zpt(ji,jj,jl) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1)
END_2D
END DO
DO jl = 1, jpl !-- flux in y-direction
DO_2D( nn_hls-1, nn_hls-1, nn_hls, nn_hls-1 )
pfv_ho(ji,jj,jl) = 0.5_wp * pv(ji,jj) * ( zpt(ji,jj,jl) + zpt(ji,jj+1,jl) )
END_2D
END DO
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_y( pdt, pv, pt, pfv_ups, pfv_ho )
ELSE !== even ice time step: adv_y then adv_x ==!
!
DO jl = 1, jpl !-- flux in y-direction
DO_2D( nn_hls, nn_hls, nn_hls, nn_hls-1 )
pfv_ho(ji,jj,jl) = 0.5_wp * pv(ji,jj) * ( pt(ji,jj,jl) + pt(ji,jj+1,jl) )
END_2D
END DO
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_y( pdt, pv, pt, pfv_ups, pfv_ho )
!
DO jl = 1, jpl !-- first guess of tracer from v-flux
DO_2D( nn_hls, nn_hls, nn_hls-1, nn_hls-1 )
ztra = - ( pfv_ho(ji,jj,jl) - pfv_ho(ji,jj-1,jl) ) &
& + ( pv (ji,jj ) - pv (ji,jj-1 ) ) * pt(ji,jj,jl) * (1.-pamsk)
!
zpt(ji,jj,jl) = ( pt(ji,jj,jl) + ztra * pdt * r1_e1e2t(ji,jj) ) * tmask(ji,jj,1)
END_2D
END DO
!
DO jl = 1, jpl !-- flux in x-direction
DO_2D( nn_hls, nn_hls-1, nn_hls-1, nn_hls-1 )
pfu_ho(ji,jj,jl) = 0.5_wp * pu(ji,jj) * ( zpt(ji,jj,jl) + zpt(ji+1,jj,jl) )
END_2D
END DO
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_x( pdt, pu, pt, pfu_ups, pfu_ho )
ENDIF
IF( np_limiter == 1 ) CALL nonosc_ice( pamsk, pdt, pu, pv, pt, pt_ups, pfu_ups, pfv_ups, pfu_ho, pfv_ho )
ENDIF
END SUBROUTINE cen2
SUBROUTINE macho( pamsk, kn_umx, jt, kt, pdt, pt, pu, pv, pubox, pvbox, pt_ups, pfu_ups, pfv_ups, pfu_ho, pfv_ho )
!!---------------------------------------------------------------------
!! *** ROUTINE macho ***
!!
!! ** Purpose : compute the high order fluxes using Ultimate-Macho scheme
!!
!! ** Method : ...
!!
!! Reference : Leonard, B.P., 1991, Comput. Methods Appl. Mech. Eng., 88, 17-74.
!!----------------------------------------------------------------------
REAL(wp) , INTENT(in ) :: pamsk ! advection of concentration (1) or other tracers (0)
INTEGER , INTENT(in ) :: kn_umx ! order of the scheme (1-5=UM or 20=CEN2)
INTEGER , INTENT(in ) :: jt ! number of sub-iteration
INTEGER , INTENT(in ) :: kt ! number of iteration
REAL(wp) , INTENT(in ) :: pdt ! tracer time-step
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pt ! tracer fields
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pu, pv ! 2 ice velocity components
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pubox, pvbox ! upstream velocity
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pt_ups ! upstream guess of tracer
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pfu_ups, pfv_ups ! upstream fluxes
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out) :: pfu_ho, pfv_ho ! high order fluxes
!
INTEGER :: ji, jj, jl ! dummy loop indices
REAL(wp), DIMENSION(jpi,jpj,jpl) :: zt_u, zt_v, zpt
!!----------------------------------------------------------------------
!
IF( MOD( (kt - 1) / nn_fsbc , 2 ) == MOD( (jt - 1) , 2 ) ) THEN !== odd ice time step: adv_x then adv_y ==!
!
! !-- ultimate interpolation of pt at u-point --!
CALL ultimate_x( nn_hls, pamsk, kn_umx, pdt, pt, pu, zt_u, pfu_ho )
! !-- limiter in x --!
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_x( pdt, pu, pt, pfu_ups, pfu_ho )
! !-- advective form update in zpt --!
DO jl = 1, jpl
DO_2D( 0, 0, nn_hls, nn_hls )
zpt(ji,jj,jl) = ( pt(ji,jj,jl) - ( pubox(ji,jj ) * ( zt_u(ji,jj,jl) - zt_u(ji-1,jj,jl) ) * r1_e1t (ji,jj) &
& + pt (ji,jj,jl) * ( pu (ji,jj ) - pu (ji-1,jj ) ) * r1_e1e2t(ji,jj) &
& * pamsk &
& ) * pdt ) * tmask(ji,jj,1)
END_2D
END DO
!
! !-- ultimate interpolation of pt at v-point --!
IF( ll_hoxy ) THEN
CALL ultimate_y( 0, pamsk, kn_umx, pdt, zpt, pv, zt_v, pfv_ho )
ELSE
CALL ultimate_y( 0, pamsk, kn_umx, pdt, pt , pv, zt_v, pfv_ho )
ENDIF
! !-- limiter in y --!
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_y( pdt, pv, pt, pfv_ups, pfv_ho )
!
!
ELSE !== even ice time step: adv_y then adv_x ==!
!
! !-- ultimate interpolation of pt at v-point --!
CALL ultimate_y( nn_hls, pamsk, kn_umx, pdt, pt, pv, zt_v, pfv_ho )
! !-- limiter in y --!
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_y( pdt, pv, pt, pfv_ups, pfv_ho )
! !-- advective form update in zpt --!
DO jl = 1, jpl
DO_2D( nn_hls, nn_hls, 0, 0 )
zpt(ji,jj,jl) = ( pt(ji,jj,jl) - ( pvbox(ji,jj ) * ( zt_v(ji,jj,jl) - zt_v(ji,jj-1,jl) ) * r1_e2t (ji,jj) &
& + pt (ji,jj,jl) * ( pv (ji,jj ) - pv (ji,jj-1 ) ) * r1_e1e2t(ji,jj) &
& * pamsk &
& ) * pdt ) * tmask(ji,jj,1)
END_2D
END DO
!
! !-- ultimate interpolation of pt at u-point --!
IF( ll_hoxy ) THEN
CALL ultimate_x( 0, pamsk, kn_umx, pdt, zpt, pu, zt_u, pfu_ho )
ELSE
CALL ultimate_x( 0, pamsk, kn_umx, pdt, pt , pu, zt_u, pfu_ho )
ENDIF
! !-- limiter in x --!
IF( np_limiter == 2 .OR. np_limiter == 3 ) CALL limiter_x( pdt, pu, pt, pfu_ups, pfu_ho )
!
ENDIF
IF( np_limiter == 1 ) CALL nonosc_ice( pamsk, pdt, pu, pv, pt, pt_ups, pfu_ups, pfv_ups, pfu_ho, pfv_ho )
!
END SUBROUTINE macho
SUBROUTINE ultimate_x( kloop, pamsk, kn_umx, pdt, pt, pu, pt_u, pfu_ho )
!!---------------------------------------------------------------------
!! *** ROUTINE ultimate_x ***
!!
!! ** Purpose : compute tracer at u-points
!!
!! ** Method : ...
!!
!! Reference : Leonard, B.P., 1991, Comput. Methods Appl. Mech. Eng., 88, 17-74.
!!----------------------------------------------------------------------
INTEGER , INTENT(in ) :: kloop ! either 0 or nn_hls depending on the order of the call
REAL(wp) , INTENT(in ) :: pamsk ! advection of concentration (1) or other tracers (0)
INTEGER , INTENT(in ) :: kn_umx ! order of the scheme (1-5=UM or 20=CEN2)
REAL(wp) , INTENT(in ) :: pdt ! tracer time-step
REAL(wp), DIMENSION(:,: ) , INTENT(in ) :: pu ! ice i-velocity component
REAL(wp), DIMENSION(:,:,:) , INTENT(in ) :: pt ! tracer fields
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out) :: pt_u ! tracer at u-point
REAL(wp), DIMENSION(jpi,jpj,jpl), INTENT( out) :: pfu_ho ! high order flux
!
INTEGER :: ji, jj, jl ! dummy loop indices
REAL(wp) :: zcu, zdx2, zdx4 ! - -
REAL(wp), DIMENSION(jpi,jpj,jpl) :: ztu1, ztu2, ztu3, ztu4
!!----------------------------------------------------------------------
!
! !-- Laplacian in i-direction --!
DO jl = 1, jpl
DO_2D( nn_hls, nn_hls-1, kloop, kloop ) ! First derivative (gradient)
ztu1(ji,jj,jl) = ( pt(ji+1,jj,jl) - pt(ji,jj,jl) ) * r1_e1u(ji,jj) * umask(ji,jj,1)
END_2D
DO_2D( nn_hls-1, nn_hls-1, kloop, kloop ) ! Second derivative (Laplacian)
ztu2(ji,jj,jl) = ( ztu1(ji,jj,jl) - ztu1(ji-1,jj,jl) ) * r1_e1t(ji,jj)
END_2D
!!$ DO jj = 2, jpjm1 ! First derivative (gradient)
!!$ DO ji = 1, jpim1
!!$ ztu1(ji,jj,jl) = ( pt(ji+1,jj,jl) - pt(ji,jj,jl) ) * r1_e1u(ji,jj) * umask(ji,jj,1)
!!$ END DO
!!$ ! ! Second derivative (Laplacian)
!!$ DO ji = 2, jpim1
!!$ ztu2(ji,jj,jl) = ( ztu1(ji,jj,jl) - ztu1(ji-1,jj,jl) ) * r1_e1t(ji,jj)
!!$ END DO
!!$ END DO
END DO
IF( nn_hls == 1 ) CALL lbc_lnk( 'icedyn_adv_umx', ztu2, 'T', 1.0_wp )
!
! !-- BiLaplacian in i-direction --!
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop ) ! Third derivative
ztu3(ji,jj,jl) = ( ztu2(ji+1,jj,jl) - ztu2(ji,jj,jl) ) * r1_e1u(ji,jj) * umask(ji,jj,1)
END_2D
DO_2D( 0, 0, kloop, kloop ) ! Fourth derivative
ztu4(ji,jj,jl) = ( ztu3(ji,jj,jl) - ztu3(ji-1,jj,jl) ) * r1_e1t(ji,jj)
END_2D
!!$ DO jj = 2, jpjm1 ! Third derivative
!!$ DO ji = 1, jpim1
!!$ ztu3(ji,jj,jl) = ( ztu2(ji+1,jj,jl) - ztu2(ji,jj,jl) ) * r1_e1u(ji,jj) * umask(ji,jj,1)
!!$ END DO
!!$ ! ! Fourth derivative
!!$ DO ji = 2, jpim1
!!$ ztu4(ji,jj,jl) = ( ztu3(ji,jj,jl) - ztu3(ji-1,jj,jl) ) * r1_e1t(ji,jj)
!!$ END DO
!!$ END DO
END DO
!
!
SELECT CASE (kn_umx )
!
CASE( 1 ) !== 1st order central TIM ==! (Eq. 21)
!
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop )
pt_u(ji,jj,jl) = 0.5_wp * umask(ji,jj,1) * ( pt(ji+1,jj,jl) + pt(ji,jj,jl) &
& - SIGN( 1._wp, pu(ji,jj) ) * ( pt(ji+1,jj,jl) - pt(ji,jj,jl) ) )
END_2D
END DO
!
CASE( 2 ) !== 2nd order central TIM ==! (Eq. 23)
!
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop )
zcu = pu(ji,jj) * r1_e2u(ji,jj) * pdt * r1_e1u(ji,jj)
pt_u(ji,jj,jl) = 0.5_wp * umask(ji,jj,1) * ( pt(ji+1,jj,jl) + pt(ji,jj,jl) &
& - zcu * ( pt(ji+1,jj,jl) - pt(ji,jj,jl) ) )
END_2D
END DO
!
CASE( 3 ) !== 3rd order central TIM ==! (Eq. 24)
!
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop )
zcu = pu(ji,jj) * r1_e2u(ji,jj) * pdt * r1_e1u(ji,jj)
zdx2 = e1u(ji,jj) * e1u(ji,jj)
!!rachid zdx2 = e1u(ji,jj) * e1t(ji,jj)
pt_u(ji,jj,jl) = 0.5_wp * umask(ji,jj,1) * ( ( pt (ji+1,jj,jl) + pt (ji,jj,jl) &
& - zcu * ( pt (ji+1,jj,jl) - pt (ji,jj,jl) ) ) &
& + z1_6 * zdx2 * ( zcu*zcu - 1._wp ) * ( ztu2(ji+1,jj,jl) + ztu2(ji,jj,jl) &
& - SIGN( 1._wp, zcu ) * ( ztu2(ji+1,jj,jl) - ztu2(ji,jj,jl) ) ) )
END_2D
END DO
!
CASE( 4 ) !== 4th order central TIM ==! (Eq. 27)
!
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop )
zcu = pu(ji,jj) * r1_e2u(ji,jj) * pdt * r1_e1u(ji,jj)
zdx2 = e1u(ji,jj) * e1u(ji,jj)
!!rachid zdx2 = e1u(ji,jj) * e1t(ji,jj)
pt_u(ji,jj,jl) = 0.5_wp * umask(ji,jj,1) * ( ( pt (ji+1,jj,jl) + pt (ji,jj,jl) &
& - zcu * ( pt (ji+1,jj,jl) - pt (ji,jj,jl) ) ) &
& + z1_6 * zdx2 * ( zcu*zcu - 1._wp ) * ( ztu2(ji+1,jj,jl) + ztu2(ji,jj,jl) &
& - 0.5_wp * zcu * ( ztu2(ji+1,jj,jl) - ztu2(ji,jj,jl) ) ) )
END_2D
END DO
!
CASE( 5 ) !== 5th order central TIM ==! (Eq. 29)
!
CALL lbc_lnk( 'icedyn_adv_umx', ztu4, 'T', 1.0_wp )
!
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop )
zcu = pu(ji,jj) * r1_e2u(ji,jj) * pdt * r1_e1u(ji,jj)
zdx2 = e1u(ji,jj) * e1u(ji,jj)
!!rachid zdx2 = e1u(ji,jj) * e1t(ji,jj)
zdx4 = zdx2 * zdx2
pt_u(ji,jj,jl) = 0.5_wp * umask(ji,jj,1) * ( ( pt (ji+1,jj,jl) + pt (ji,jj,jl) &
& - zcu * ( pt (ji+1,jj,jl) - pt (ji,jj,jl) ) ) &
& + z1_6 * zdx2 * ( zcu*zcu - 1._wp ) * ( ztu2(ji+1,jj,jl) + ztu2(ji,jj,jl) &
& - 0.5_wp * zcu * ( ztu2(ji+1,jj,jl) - ztu2(ji,jj,jl) ) ) &
& + z1_120 * zdx4 * ( zcu*zcu - 1._wp ) * ( zcu*zcu - 4._wp ) * ( ztu4(ji+1,jj,jl) + ztu4(ji,jj,jl) &
& - SIGN( 1._wp, zcu ) * ( ztu4(ji+1,jj,jl) - ztu4(ji,jj,jl) ) ) )
END_2D
END DO
!
END SELECT
!
! if pt at u-point is negative then use the upstream value
! this should not be necessary if a proper sea-ice mask is set in Ultimate
! to degrade the order of the scheme when necessary (for ex. at the ice edge)
IF( ll_neg ) THEN
DO jl = 1, jpl
DO_2D( 1, 0, kloop, kloop )
IF( pt_u(ji,jj,jl) < 0._wp .OR. ( imsk_small(ji,jj,jl) == 0 .AND. pamsk == 0. ) ) THEN
pt_u(ji,jj,jl) = 0.5_wp * umask(ji,jj,1) * ( pt(ji+1,jj,jl) + pt(ji,jj,jl) &
& - SIGN( 1._wp, pu(ji,jj) ) * ( pt(ji+1,jj,jl) - pt(ji,jj,jl) ) )
ENDIF
END_2D
END DO
ENDIF
! !-- High order flux in i-direction --!
DO jl = 1, jpl
DO_2D( 1, 0, 0, 0 )
pfu_ho(ji,jj,jl) = pu(ji,jj) * pt_u(ji,jj,jl)
END_2D
END DO
!
END SUBROUTINE ultimate_x