-
Notifications
You must be signed in to change notification settings - Fork 69
Expand file tree
/
Copy pathicb_dyn.F90
More file actions
1225 lines (995 loc) · 41 KB
/
icb_dyn.F90
File metadata and controls
1225 lines (995 loc) · 41 KB
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 iceberg_dynamics
USE MOD_MESH
use MOD_PARTIT
USE MOD_PARSUP
use MOD_ICE
USE MOD_DYN
use iceberg_params
use iceberg_element
!use iceberg_step
implicit none
public :: iceberg_dyn
public :: iceberg_frozen
public :: iceberg_acceleration
public :: compute_areas
public :: iceberg_average_andkeel
public :: iceberg_avvelo
contains
!==============================================================================
! calculates basically the new iceberg velocity; if melting is enabled, the
! iceberg dimensions are adjusted as well.
!
! Thomas Rackow, 29.06.2010
!==============================================================================
subroutine iceberg_dyn(mesh, partit, ice, dynamics, ib, new_u_ib, new_v_ib, u_ib, v_ib, lon,lat, depth_ib, &
height_ib, length_ib, width_ib, iceberg_elem, &
mass_ib, Ci, Ca, Co, Cda_skin, Cdo_skin, &
rho_ice, rho_air, rho_h2o, P_sill, conc_sill, frozen_in, &
file1, file2, P_ib, conci_ib, dt_ib, lastsubstep, &
f_u_ib_old, f_v_ib_old, l_semiimplicit, &
semiimplicit_coeff, AB_coeff, file3, rho_icb)
use g_forcing_arrays !for u_wind, v_wind or u_wind_ib, v_wind_ib respectively
use o_arrays, only: Tsurf_ib, Ssurf_ib
use o_param !for dt
integer :: ib_n_lvls, m
integer, intent(IN) :: ib !current iceberg's index
real, intent(OUT) :: new_u_ib, new_v_ib
real, intent(IN) :: u_ib, v_ib
real, intent(IN) :: lon,lat !radiant
real, intent(INOUT) :: depth_ib !inout for case of melting iceberg
real, intent(INOUT) :: height_ib !inout for case of melting iceberg
real, intent(INOUT) :: length_ib !inout for case of melting iceberg
real, intent(INOUT) :: width_ib !inout for case of melting iceberg
integer, intent(IN) :: iceberg_elem !local
real, intent(IN) :: mass_ib
real, intent(IN) :: Ci, Ca, Co, Cda_skin, Cdo_skin
real, intent(IN) :: rho_ice, rho_air, rho_h2o
real, intent(IN) :: P_sill, conc_sill
real, intent(INOUT) :: frozen_in
character, intent(IN) :: file1*80, file2*80
real, intent(OUT) :: P_ib, conci_ib
real, intent(IN) :: dt_ib
logical, intent(IN) :: lastsubstep
real, intent(INOUT) :: f_u_ib_old, f_v_ib_old
logical, intent(IN) :: l_semiimplicit
real, intent(IN) :: semiimplicit_coeff, AB_coeff
!LA 2023-03-07
real, dimension(:), pointer :: hi_ib3, conci_ib3, coriolis
real, dimension(3) :: uo_keel, vo_keel, T_keel,S_keel, uo_dz, vo_dz, T_dz,S_dz!hi_ib3, conci_ib3,
real, dimension(:,:), allocatable :: arr_uo_dz, arr_vo_dz, arr_T_dz,arr_S_dz !hi_ib3, conci_ib3,
real :: uo_ib, vo_ib, ua_ib, va_ib, ui_ib, vi_ib, hi_ib, uo_skin_ib, vo_skin_ib
real, dimension(:), allocatable :: arr_uo_ib, arr_vo_ib, arr_T_ave_ib, arr_S_ave_ib
real :: Ao, Aa, Ai, Ad, fcoriolis
real :: au_ib, av_ib
real, dimension(2,2) :: SI_matrix
real, dimension(2) :: SI_velo
real :: u_ib_tmp, v_ib_tmp, normold, normnew, abs_omib, abs_omib_skin, ocean_drag
integer :: iter_ib, n, n2
real :: M_b, M_v, M_e, M_bv, sst_ib, sss_ib ! meltrates (basal, lateral, erosion, lateral 'basal'), temp. & salinity
real :: T_ave_ib, S_ave_ib, T_keel_ib, S_keel_ib
character, intent(IN) :: file3*80
real, intent(IN) :: rho_icb
type(t_ice) , intent(inout), target :: ice
type(t_mesh), intent(in) , target :: mesh
type(t_partit), intent(inout), target :: partit
type(t_dyn) , intent(inout), target :: dynamics
#include "associate_part_def.h"
#include "associate_mesh_def.h"
#include "associate_part_ass.h"
#include "associate_mesh_ass.h"
n2=elem2D_nodes(1,iceberg_elem)
allocate(arr_uo_dz(3,nlevels_nod2D(n2)))
allocate(arr_vo_dz(3,nlevels_nod2D(n2)))
allocate(arr_T_dz(3,nlevels_nod2D(n2)))
allocate(arr_S_dz(3,nlevels_nod2D(n2)))
arr_uo_dz = 0.0
arr_vo_dz = 0.0
arr_T_dz = 0.0
arr_S_dz = 0.0
allocate(arr_uo_ib(nlevels_nod2D(n2)))
allocate(arr_vo_ib(nlevels_nod2D(n2)))
allocate(arr_T_ave_ib(nlevels_nod2D(n2)))
allocate(arr_S_ave_ib(nlevels_nod2D(n2)))
arr_uo_ib = 0.0
arr_vo_ib = 0.0
arr_T_ave_ib = 0.0
arr_S_ave_ib = 0.0
!OCEAN VELOCITIES:
! - (uo_ib, vo_ib) : integrated mean velocity at location of iceberg
! - (uo_skin_ib, vo_skin_ib) : velocity below the draft of the iceberg
! call iceberg_avvelo_ufkeel(uo_dz,vo_dz, uo_keel,vo_keel, depth_ib,iceberg_elem)
call iceberg_average_andkeel(mesh, partit, dynamics, uo_dz,vo_dz, uo_keel,vo_keel, T_dz,S_dz, T_keel,S_keel, depth_ib,iceberg_elem, ib)
!OCEANIC VELOCITY uo_ib, vo_ib
call FEM_3eval(mesh, partit, uo_ib,vo_ib,lon,lat,uo_dz,vo_dz,iceberg_elem)
call iceberg_levelwise_andkeel(mesh, partit, dynamics, arr_uo_dz,arr_vo_dz, uo_keel,vo_keel, arr_T_dz,arr_S_dz, T_keel,S_keel, depth_ib,iceberg_elem, ib, ib_n_lvls)
do n=1, ib_n_lvls
call FEM_3eval(mesh, partit, arr_uo_ib(n),arr_vo_ib(n),lon,lat,arr_uo_dz(:,n),arr_vo_dz(:,n),iceberg_elem)
call FEM_3eval(mesh, partit, arr_T_ave_ib(n),arr_S_ave_ib(n),lon,lat,arr_T_dz(:,n),arr_S_dz(:,n),iceberg_elem)
end do
call FEM_3eval(mesh, partit, uo_skin_ib,vo_skin_ib,lon,lat,uo_keel,vo_keel,iceberg_elem)
!TEMPERATURE AND SALINITY:
! - T_ave_ib, S_ave_ib : Mean T & S (integrated) at location of iceberg
! - T_keel_ib, S_keel_ib : T & S below the draft of the iceberg (depth_ib)
call FEM_3eval(mesh, partit, T_keel_ib,S_keel_ib,lon,lat,T_keel,S_keel,iceberg_elem)
!ATMOSPHERIC VELOCITY ua_ib, va_ib
call FEM_eval(mesh, partit, ua_ib,va_ib,lon,lat,u_wind_ib,v_wind_ib,iceberg_elem)
!ICE VELOCITY ui_ib, vi_ib
call FEM_eval(mesh, partit, ui_ib,vi_ib,lon,lat,ice%uice_ib,ice%vice_ib,iceberg_elem)
!ICE THICKNESS (CONCENTRATION) hi_ib, conci_ib
hi_ib3 => ice%data(size(ice%data))%values(:) !ice%m_ice_ib(tmp_arr)
conci_ib3 => ice%data(size(ice%data)-1)%values(:) !ice%a_ice_ib(tmp_arr)
call FEM_3eval(mesh, partit, hi_ib,conci_ib,lon,lat,hi_ib3,conci_ib3,iceberg_elem)
P_ib = 20000. * hi_ib * exp(-20.*(1-conci_ib))
call compute_areas(Ao, Aa, Ai, Ad, depth_ib, &
height_ib, length_ib, width_ib, hi_ib)
coriolis => mesh%coriolis(:)
fcoriolis = coriolis(iceberg_elem) * coriolis_scale(ib)
call iceberg_acceleration( mesh, partit, dynamics, ib, au_ib, av_ib, Ao, Aa, Ai, Ad, &
uo_ib,vo_ib, uo_skin_ib, vo_skin_ib, &
ua_ib,va_ib, ui_ib,vi_ib, &
u_ib, v_ib, mass_ib, fcoriolis, &
Ci, Ca, Co, Cda_skin, Cdo_skin, &
rho_ice, rho_air, rho_h2o, length_ib, &
iceberg_elem, conci_ib, file1, file2, &
lon, lat, lastsubstep, f_u_ib_old, &
f_v_ib_old, l_semiimplicit, &
semiimplicit_coeff, AB_coeff )
!========================THERMODYNAMICS============================
if(l_melt) then
call FEM_eval(mesh, partit, sst_ib,sss_ib,lon,lat,Tsurf_ib,Ssurf_ib,iceberg_elem)
call iceberg_meltrates(partit, mesh, M_b, M_v, M_e, M_bv, &
u_ib,v_ib, arr_uo_ib,arr_vo_ib, ua_ib,va_ib, &
sst_ib, length_ib, conci_ib, &
uo_skin_ib, vo_skin_ib, T_keel_ib, S_keel_ib, depth_ib, height_ib, &
arr_T_ave_ib, arr_S_ave_ib, ib, rho_icb, uo_ib, vo_ib, ib_n_lvls, iceberg_elem, nlevels_nod2D(n2))
call iceberg_newdimensions(partit, ib, depth_ib,height_ib,length_ib,width_ib,M_b,M_v,M_e,M_bv, &
rho_h2o, rho_icb, file3)
end if
!====================END OF THERMODYNAMICS=========================
new_u_ib = u_ib + au_ib * dt_ib
new_v_ib = v_ib + av_ib * dt_ib
if (l_semiimplicit) then !a matrix multiplication is to be performed
!for semiimpl. coriolis term and implicit
!water drag
abs_omib = sqrt( (uo_ib - u_ib)**2 + (vo_ib - v_ib)**2 )
abs_omib_skin = sqrt( (uo_skin_ib - u_ib)**2 + (vo_skin_ib - v_ib)**2 )
ocean_drag = (0.5 * Co * rho_h2o * Ao * abs_omib + rho_h2o * Cdo_skin * Ad &
* abs_omib_skin)/mass_ib
SI_matrix(1,1) = 1. + dt_ib*ocean_drag
SI_matrix(1,2) = dt_ib*fcoriolis*semiimplicit_coeff
SI_matrix(2,1) =-SI_matrix(1,2)
SI_matrix(2,2) = SI_matrix(1,1)
SI_matrix = (1./( SI_matrix(2,2)**2 + SI_matrix(1,2)**2 )) * SI_matrix
!new velocity
SI_velo = MATMUL(SI_matrix, (/ new_u_ib, new_v_ib /))
!for P_sill calculate "effective acceleration" in case of semi-impl. scheme by
!(new velo - old velo)/dt for au_ib (av_ib) isn't the real acceleration
!au_ib = (SI_velo(1) - u_ib)/dt_ib
!av_ib = (SI_velo(2) - v_ib)/dt_ib
!now the velocity can be updated
new_u_ib = SI_velo(1)
new_v_ib = SI_velo(2)
else !compute only water drag implicitly, coriolis: AB
abs_omib = sqrt( (uo_ib - u_ib)**2 + (vo_ib - v_ib)**2 )
abs_omib_skin = sqrt( (uo_skin_ib - u_ib)**2 + (vo_skin_ib - v_ib)**2 )
ocean_drag = (0.5 * Co * rho_h2o * Ao * abs_omib + rho_h2o * Cdo_skin * Ad &
* abs_omib_skin)/mass_ib
SI_matrix(1,1) = 1. + dt_ib*ocean_drag
SI_matrix(1,2) = 0.0
SI_matrix(2,1) = 0.0
SI_matrix(2,2) = SI_matrix(1,1)
SI_matrix = (1./(SI_matrix(2,2)**2)) * SI_matrix
!new velocity
SI_velo = MATMUL(SI_matrix, (/ new_u_ib, new_v_ib /))
!now the velocity can be updated
new_u_ib = SI_velo(1)
new_v_ib = SI_velo(2)
end if !matrix-multiplication
!icebergs may be frozen in: .true. or .false.
call iceberg_frozen(frozen_in, P_sill, P_ib, conc_sill, conci_ib, ib)
new_u_ib = (1-frozen_in) * new_u_ib + frozen_in * ui_ib
new_v_ib = (1-frozen_in) * new_v_ib + frozen_in * vi_ib
end subroutine iceberg_dyn
!***************************************************************************************************************************
!***************************************************************************************************************************
subroutine iceberg_frozen(festgefroren, P_sill, P_ib, conc_sill, conci_ib, ib)
!use iceberg_params, only : l_freeze !use 'capturing mechanism' of sea ice?
implicit none
real, intent(OUT) :: festgefroren
real, intent(IN) :: P_sill, P_ib
real, intent(IN) :: conc_sill, conci_ib
integer, intent(IN) :: ib
festgefroren=0.0
if( l_freeze(ib) ) then
festgefroren= factor(P_ib,P_sill,P_sill-3000.) &
* factor(conci_ib,conc_sill,conc_sill-0.04)
end if
contains
!=================================================================
! ... a function for linear transition zone between "frozen in"
! and "not frozen in", where sill & zone are defined as below
!
! sill
! |
! /--------- 1
! 0 ________/
! |
! zone and value is P_ib or conci_ib.
!
!=================================================================
real function factor(value, sill,zone)
implicit none
real, intent(IN) :: value
real, intent(IN) :: sill
real, intent(IN) :: zone
if(value < zone) then
factor=0.
else if(value > sill) then
factor=1.
else
factor=(value-zone)/(sill-zone)
end if
end function factor
end subroutine iceberg_frozen
!***************************************************************************************************************************
!***************************************************************************************************************************
subroutine iceberg_acceleration(mesh, partit, dynamics, ib, au_ib, av_ib, Ao, Aa, Ai, Ad, &
uo_ib,vo_ib, uo_skin_ib, vo_skin_ib, &
ua_ib,va_ib, ui_ib,vi_ib, &
u_ib, v_ib, mass_ib, fcoriolis, &
Ci, Ca, Co, Cda_skin, Cdo_skin, &
rho_ice, rho_air, rho_h2o, length_ib, &
iceberg_elem, conci_ib, file1, file2, &
lon_rad, lat_rad, output, f_u_ib_old, &
f_v_ib_old, l_semiimplicit, &
semiimplicit_coeff, AB_coeff )
use o_param !for g
use g_config !for istep
use MOD_PARTIT !for mype
use g_rotate_grid, only: vector_r2g, vector_g2r
!use iceberg_params, only: l_wave, l_tides, l_geo_out, surfslop_scale, ascii_out
use MOD_MESH
use MOD_DYN
implicit none
integer, intent(IN) :: ib
real, intent(OUT) :: au_ib, av_ib
real, intent(IN) :: Ao, Aa, Ai, Ad
real, intent(IN) :: uo_ib,vo_ib, uo_skin_ib, vo_skin_ib, ua_ib,va_ib, ui_ib,vi_ib, u_ib,v_ib
real, intent(IN) :: mass_ib, fcoriolis
real, intent(IN) :: Ci, Ca, Co, Cda_skin, Cdo_skin
real, intent(IN) :: rho_ice, rho_air, rho_h2o, length_ib
integer, intent(IN) :: iceberg_elem
real, intent(IN) :: conci_ib
character, intent(IN) :: file1*80, file2*80
real, intent(IN) :: lon_rad, lat_rad
logical, intent(IN) :: output
real, intent(INOUT):: f_u_ib_old, f_v_ib_old
logical, intent(IN) :: l_semiimplicit
real, intent(IN) :: semiimplicit_coeff, AB_coeff
real :: vel_atm, wave_amplitude, direction_u, direction_v
real :: abs_omib, abs_amib, abs_imib, abs_omib_skin
real :: ocean_drag_u, ocean_skin_u, air_drag_u, air_skin_u
real :: ice_drag_u, wave_radiation_u, surface_slope_u, slope_tides_u
real :: ocean_drag_v, ocean_skin_v, air_drag_v, air_skin_v
real :: ice_drag_v, wave_radiation_v, surface_slope_v, slope_tides_v
real, dimension(2) :: nablaeta
real, dimension(8) :: accs_u_out, accs_v_out
real, dimension(4) :: vels_u_out, vels_v_out
real :: oneminus_AB !, test1, test2
integer :: i, istep
type(t_mesh), intent(in) , target :: mesh
type(t_partit), intent(inout), target :: partit
type(t_dyn) , intent(inout), target :: dynamics
#include "associate_part_def.h"
#include "associate_mesh_def.h"
#include "associate_part_ass.h"
#include "associate_mesh_ass.h"
!estimate wave height at the icebergs location (Bigg et al., 1997),
!so wave_amplitude = 0.5 * wave_height = 0.5 * const. * abs(atm velo)**2
vel_atm = sqrt(ua_ib**2 + va_ib**2)
wave_amplitude = 0.5 * 0.02025 * vel_atm**2
!assume that waves have same direction as the winds
direction_u = ua_ib / vel_atm
direction_v = va_ib / vel_atm
!absolute values of relative velocities
abs_omib = sqrt( (uo_ib - u_ib)**2 + (vo_ib - v_ib)**2 )
abs_amib = sqrt( (ua_ib - u_ib)**2 + (va_ib - v_ib)**2 )
abs_imib = sqrt( (ui_ib - u_ib)**2 + (vi_ib - v_ib)**2 )
abs_omib_skin = sqrt( (uo_skin_ib - u_ib)**2 + (vo_skin_ib - v_ib)**2 )
! u-components
ocean_drag_u = (0.5 * Co * rho_h2o * Ao &
* abs_omib * uo_ib)/mass_ib !calculate part of it implicitly
ocean_skin_u = (rho_h2o * Cdo_skin * Ad &
* abs_omib_skin * uo_skin_ib)/mass_ib !calculate part of it implicitly
air_drag_u = (0.5 * Ca * rho_air * Aa &
* abs_amib * (ua_ib - u_ib))/mass_ib
air_skin_u = (rho_air * Cda_skin * Ad &
* abs_amib * (ua_ib - u_ib))/mass_ib
ice_drag_u = (0.5 * Ci * rho_ice * Ai &
* abs_imib * (ui_ib - u_ib))/mass_ib
if(l_wave(ib)) then
wave_radiation_u = 1./4. * rho_h2o * g * length_ib &
* wave_amplitude**2 * direction_u /mass_ib
else
wave_radiation_u = 0.0
end if
!use gradient smoothing for surface slope term
call mean_gradient(mesh, partit, dynamics, iceberg_elem, lon_rad, lat_rad, nablaeta)
!additional surface slope due to tides
if(l_tides) then
slope_tides_u = 0.0 !-g * sum( opbnd_z_tide(elem2D_nodes(:,iceberg_elem)) * bafux_2D(:,iceberg_elem) )
else
slope_tides_u = 0.0
end if
surface_slope_u = (-g * nablaeta(1) + slope_tides_u) * surfslop_scale(ib) !default scaling is 1.0
! v-components
ocean_drag_v = (0.5 * Co * rho_h2o * Ao &
* abs_omib * vo_ib)/mass_ib !calculate part of it implicitly
ocean_skin_v = (rho_h2o * Cdo_skin * Ad &
* abs_omib_skin * vo_skin_ib)/mass_ib !calculate part of it implicitly
air_drag_v = (0.5 * Ca * rho_air * Aa &
* abs_amib * (va_ib - v_ib))/mass_ib
air_skin_v = (rho_air * Cda_skin * Ad &
* abs_amib * (va_ib - v_ib))/mass_ib
ice_drag_v = (0.5 * Ci * rho_ice * Ai &
* abs_imib * (vi_ib - v_ib))/mass_ib
if(l_wave(ib)) then
wave_radiation_v = 1./4. * rho_h2o * g * length_ib &
* wave_amplitude**2 * direction_v /mass_ib
else
wave_radiation_v = 0.0
end if
!additional surface slope due to tides
if(l_tides) then
slope_tides_v = 0.0 !-g * sum( opbnd_z_tide(elem2D_nodes(:,iceberg_elem)) * bafuy_2D(:,iceberg_elem) )
else
slope_tides_v = 0.0
end if
surface_slope_v = (-g * nablaeta(2) + slope_tides_v) * surfslop_scale(ib) !default scaling is 1.0
if (l_semiimplicit) then !USE (SEMI-)IMPLICIT SCHEME for coriolis term
if (conci_ib .GT. 0.15) then
au_ib = ocean_drag_u &
+ ocean_skin_u &
+ air_drag_u &
+ air_skin_u &
+ ice_drag_u &
+ wave_radiation_u &
+ surface_slope_u &
+ (1.-semiimplicit_coeff)*fcoriolis*v_ib
av_ib = ocean_drag_v &
+ ocean_skin_v &
+ air_drag_v &
+ air_skin_v &
+ ice_drag_v &
+ wave_radiation_v &
+ surface_slope_v &
- (1.-semiimplicit_coeff)*fcoriolis*u_ib
else
au_ib = ocean_drag_u &
+ ocean_skin_u &
+ air_drag_u &
+ air_skin_u &
+ wave_radiation_u &
+ surface_slope_u &
+ (1.-semiimplicit_coeff)*fcoriolis*v_ib
av_ib = ocean_drag_v &
+ ocean_skin_v &
+ air_drag_v &
+ air_skin_v &
+ wave_radiation_v &
+ surface_slope_v &
- (1.-semiimplicit_coeff)*fcoriolis*u_ib
end if
else !USE ADAMS-BASHFORTH SCHEME for coriolis
oneminus_AB= 1. - AB_coeff
if (conci_ib .GT. 0.15) then
au_ib = ocean_drag_u &
+ ocean_skin_u &
+ air_drag_u &
+ air_skin_u &
+ ice_drag_u &
+ wave_radiation_u &
+ surface_slope_u &
+ AB_coeff*fcoriolis*v_ib+oneminus_AB*f_v_ib_old
av_ib = ocean_drag_v &
+ ocean_skin_v &
+ air_drag_v &
+ air_skin_v &
+ ice_drag_v &
+ wave_radiation_v &
+ surface_slope_v &
- (AB_coeff*fcoriolis*u_ib) - (oneminus_AB*f_u_ib_old)
else
au_ib = ocean_drag_u &
+ ocean_skin_u &
+ air_drag_u &
+ air_skin_u &
+ wave_radiation_u &
+ surface_slope_u &
+ AB_coeff*fcoriolis*v_ib+oneminus_AB*f_v_ib_old
av_ib = ocean_drag_v &
+ ocean_skin_v &
+ air_drag_v &
+ air_skin_v &
+ wave_radiation_v &
+ surface_slope_v &
- (AB_coeff*fcoriolis*u_ib) - (oneminus_AB*f_u_ib_old)
end if
!save f*velocity for A.B. scheme before it is updated
f_u_ib_old=fcoriolis*u_ib
f_v_ib_old=fcoriolis*v_ib
end if !use semiimplicit scheme or AB method?
! !if(.false.) then
! if(output .AND. ascii_out) then
!
! accs_u_out(1) = ocean_drag_u - (0.5 * Co * rho_h2o * Ao * abs_omib * u_ib)/mass_ib
! accs_u_out(2) = ocean_skin_u - (rho_h2o * Cdo_skin * Ad * abs_omib_skin * u_ib)/mass_ib
! accs_u_out(3) = air_drag_u
! accs_u_out(4) = air_skin_u
! accs_u_out(5) = ice_drag_u
! accs_u_out(6) = wave_radiation_u
! accs_u_out(7) = surface_slope_u
! accs_u_out(8) = fcoriolis*v_ib
! vels_u_out(1) = uo_ib
! vels_u_out(2) = uo_skin_ib
! vels_u_out(3) = ua_ib
! vels_u_out(4) = ui_ib
!
! accs_v_out(1) = ocean_drag_v - (0.5 * Co * rho_h2o * Ao * abs_omib * v_ib)/mass_ib
! accs_v_out(2) = ocean_skin_v - (rho_h2o * Cdo_skin * Ad * abs_omib_skin * v_ib)/mass_ib
! accs_v_out(3) = air_drag_v
! accs_v_out(4) = air_skin_v
! accs_v_out(5) = ice_drag_v
! accs_v_out(6) = wave_radiation_v
! accs_v_out(7) = surface_slope_v
! accs_v_out(8) = -fcoriolis*u_ib
! vels_v_out(1) = vo_ib
! vels_v_out(2) = vo_skin_ib
! vels_v_out(3) = va_ib
! vels_v_out(4) = vi_ib
!
! if(l_geo_out) then
! do i=1,8
! call vector_r2g(accs_u_out(i), accs_v_out(i), lon_rad, lat_rad, 0)
! end do
! do i=1,4
! call vector_r2g(vels_u_out(i), vels_v_out(i), lon_rad, lat_rad, 0)
! end do
! end if
!
! !open(unit=icbID,file=file1,position='append')
! !write(icbID,'(2I,12e15.7)') mype, istep, &
! ! accs_u_out(1), accs_u_out(2), &
! ! accs_u_out(3), accs_u_out(4), &
! ! accs_u_out(5), accs_u_out(6), &
! ! accs_u_out(7), accs_u_out(8), &
! ! vels_u_out(1), vels_u_out(2), &
! ! vels_u_out(3), vels_u_out(4)
! !
! !
! !
! !close(icbID)
!
! !open(unit=icbID,file=file2,position='append')
! !write(icbID,'(2I,12e15.7)') mype,istep, &
! ! accs_v_out(1), accs_v_out(2), &
! ! accs_v_out(3), accs_v_out(4), &
! ! accs_v_out(5), accs_v_out(6), &
! ! accs_v_out(7), accs_v_out(8), &
! ! vels_v_out(1), vels_v_out(2), &
! ! vels_v_out(3), vels_v_out(4)
! !close(icbID)
! end if !output
end subroutine iceberg_acceleration
!***************************************************************************************************************************
!***************************************************************************************************************************
subroutine compute_areas(Ao, Aa, Ai, Ad, depth_ib, &
height_ib, length_ib, width_ib, hi_ib)
implicit none
real, intent(OUT) :: Ao, Aa, Ai, Ad
real, intent(IN) :: depth_ib, height_ib, length_ib, width_ib, hi_ib
!area of iceberg exposed to ocean, atm, seaice, horizontal area
Ao = abs(depth_ib) * length_ib
Aa = (height_ib - abs(depth_ib)) * length_ib
Ai = hi_ib * length_ib
Ad = length_ib * width_ib
end subroutine compute_areas
!***************************************************************************************************************************
!***************************************************************************************************************************
subroutine iceberg_levelwise_andkeel(mesh, partit, dynamics, uo_dz,vo_dz, uo_keel,vo_keel, T_dz,S_dz, T_keel,S_keel, depth_ib,iceberg_elem, ib, ib_n_lvls)
USE MOD_MESH
use o_param
use MOD_PARTIT
use MOD_DYN
use o_arrays, only: Tclim_ib, Sclim_ib !, UV_ib, Z_3d_n_ib
use g_clock
use g_forcing_arrays
use g_rotate_grid
implicit none
REAL, DIMENSION(:,:), allocatable, INTENT(OUT) :: uo_dz
REAL, DIMENSION(:,:), allocatable, INTENT(OUT) :: vo_dz
REAL, DIMENSION(3), INTENT(OUT) :: uo_keel
REAL, DIMENSION(3), INTENT(OUT) :: vo_keel
REAL, DIMENSION(:,:), allocatable, INTENT(OUT) :: T_dz
REAL, DIMENSION(:,:), allocatable, INTENT(OUT) :: S_dz
REAL, DIMENSION(3), INTENT(OUT) :: T_keel
REAL, DIMENSION(3), INTENT(OUT) :: S_keel
REAl, INTENT(IN) :: depth_ib
INTEGER :: ib_n_lvls_old
INTEGER, INTENT(IN) :: iceberg_elem, ib
INTEGER, INTENT(OUT) :: ib_n_lvls
INTEGER, dimension(3) :: arr_ib_n_lvls
REAL, dimension(:,:,:), pointer :: UV_ib
real :: lev_up, lev_low
integer :: m, k, n2, n_up, n_low, cavity_count, max_node_level_count
! depth over which is integrated (layer and sum)
real :: dz, ufkeel1, ufkeel2, Temkeel, Salkeel, ldepth_up, ldepth_low, dz_depth
type(t_mesh), intent(in) , target :: mesh
type(t_dyn), intent(in) , target :: dynamics
type(t_partit), intent(inout), target :: partit
#include "associate_part_def.h"
#include "associate_mesh_def.h"
#include "associate_part_ass.h"
#include "associate_mesh_ass.h"
UV_IB => dynamics%uv_ib(:,:,:)
cavity_count=0
do m=1,3
if(m==1) then
max_node_level_count = nlevels_nod2D(elem2D_nodes(m,iceberg_elem))
else
max_node_level_count = max(max_node_level_count, nlevels_nod2D(elem2D_nodes(m,iceberg_elem)))
end if
end do
allocate(uo_dz(3,max_node_level_count))
allocate(vo_dz(3,max_node_level_count))
allocate(T_dz(3,max_node_level_count))
allocate(S_dz(3,max_node_level_count))
ib_n_lvls_old = 0
ib_n_lvls = 0
arr_ib_n_lvls = 0
uo_dz = 0.0
vo_dz = 0.0
uo_keel = 0.0
vo_keel = 0.0
T_dz = 0.0
S_dz = 0.0
T_keel = 0.0
S_keel = 0.0
!LOOP: over all nodes of the iceberg element
nodeloop: do m=1, 3
!for each 2D node of the iceberg element..
n2=elem2D_nodes(m,iceberg_elem)
! LOOP: consider all neighboring pairs (n_up,n_low) of 3D nodes
! below n2..
!innerloop: do k=1, nl+1
innerloop: do k=1, nlevels_nod2D(n2)
lev_up = mesh%zbar_3d_n(k, n2)
!lev_up = mesh%Z_3d_n(k, n2)
ldepth_up = mesh%Z_3d_n(k, n2)
if( k==nlevels_nod2D(n2) ) then
lev_low = mesh%zbar_n_bot(n2)
ldepth_low = mesh%zbar_n_bot(n2)
else
lev_low = mesh%zbar_3d_n(k+1, n2)
ldepth_low = mesh%Z_3d_n(k+1, n2)
!lev_low = mesh%Z_3d_n(k+1, n2)
end if
!if( k==1 ) then
! lev_up = 0.0
!else
! lev_up = mesh%Z_3d_n_ib(k-1, n2)
! !lev_up = mesh%Z_3d_n_ib(k-1, n2)
!end if
!if( k==nlevels_nod2D(n2) ) then
! lev_low = mesh%zbar_n_bot(n2)
!else
! lev_low = mesh%Z_3d_n_ib(k, n2)
!end if
dz = abs( lev_low - lev_up )
dz_depth = abs( ldepth_low - ldepth_up )
!if( abs(lev_up)>=abs(depth_ib) ) then
! ! ...icb bottom above lev_up --> no further integration
!end if
!if( (abs(coord_nod3D(3, n_low))>abs(depth_ib)) .AND. (abs(coord_nod3D(3, n_up))>abs(depth_ib)) ) then
! write(*,*) 'INFO, k:',k,'z_up:',coord_nod3D(3, n_up),'z_lo:',coord_nod3D(3, n_low),'depth:',depth_ib,'cavity:',(cavity_flag_nod2d(elem2D_nodes(m,iceberg_elem))==1)
!end if
! if cavity node ..
if (use_cavity .AND. mesh%cavity_depth(elem2D_nodes(m,iceberg_elem)) /= 0.0 .AND. abs(depth_ib) < abs(lev_up)) then
! LA: Never go here for k=1, because abs(depth_ib)>=0.0 for all icebergs
uo_dz(m,k)=UV_ib(1,k-1,n2)*abs(depth_ib)
vo_dz(m,k)=UV_ib(2,k-1,n2)*abs(depth_ib)
uo_keel(m)=UV_ib(1,k-1,n2)
vo_keel(m)=UV_ib(2,k-1,n2)
T_dz(m,k)=Tclim_ib(k-1,n2)*abs(depth_ib)
S_dz(m,k)=Sclim_ib(k-1,n2)*abs(depth_ib)
T_keel(m)=Tclim_ib(k-1,n2)
S_keel(m)=Sclim_ib(k-1,n2) ! check those choices with RT: OK
exit innerloop
! if the lowest z coord is below the iceberg draft, exit
!else if( abs(coord_nod3D(3, n_low))>=abs(depth_ib) .AND. abs(coord_nod3D(3, n_up))<=abs(depth_ib) ) then
!****************************************************************
! LA 23.11.21 case if depth_ib<lev_up
else !# comp cav flag
if( abs(lev_low)>=abs(depth_ib) ) then !.AND. (abs(lev_up)<=abs(depth_ib)) ) then
if( abs(lev_up)<abs(depth_ib) ) then
dz = abs ( lev_up - depth_ib )
else
! LA: Never go here, when starting with k=1
dz = abs(depth_ib)
end if
!****************************************************************
if( k==1 ) then
ufkeel1 = UV_ib(1,k,n2)
ufkeel2 = UV_ib(2,k,n2)
Temkeel = Tclim_ib(k,n2)
Salkeel = Sclim_ib(k,n2)
else if( k.eq.nlevels_nod2D(n2) ) then
ufkeel1 = UV_ib(1,k-1,n2)
ufkeel2 = UV_ib(2,k-1,n2)
Temkeel = Tclim_ib(k-1,n2)
Salkeel = Sclim_ib(k-1,n2)
else
ufkeel1 = interpol1D(abs(lev_up),UV_ib(1,k-1,n2),abs(lev_low),UV_ib(1,k,n2),abs(depth_ib))
ufkeel2 = interpol1D(abs(lev_up),UV_ib(2,k-1,n2),abs(lev_low),UV_ib(2,k,n2),abs(depth_ib))
Temkeel = interpol1D(abs(lev_up),Tclim_ib(k-1,n2),abs(lev_low),Tclim_ib(k,n2),abs(depth_ib))
Salkeel = interpol1D(abs(lev_up),Sclim_ib(k-1,n2),abs(lev_low),Sclim_ib(k,n2),abs(depth_ib))
end if
uo_dz(m,k)=ufkeel1
vo_dz(m,k)=ufkeel2
T_dz(m,k)=Temkeel
S_dz(m,k)=Salkeel
uo_keel(m)=ufkeel1
vo_keel(m)=ufkeel2
T_keel(m)=Temkeel
S_keel(m)=Salkeel
arr_ib_n_lvls(m) = k
exit innerloop
!****************************************************************
! LA 23.11.21 case if lev_low==0
else if(lev_low==lev_up) then
arr_ib_n_lvls(m) = k
exit innerloop
!****************************************************************
else
if( k==1 ) then
uo_dz(m,k)=UV_ib(1,k,n2)
vo_dz(m,k)=UV_ib(2,k,n2)
T_dz(m,k)=Tclim_ib(k,n2)
S_dz(m,k)=Sclim_ib(k,n2)
elseif (k.eq.nlevels_nod2D(n2)) then ! LA 2023-08-31
! .. and sum up the layer-integrated velocities ..
! kh 08.03.21 use UV_ib buffered values here
uo_dz(m,k)=UV_ib(1,k-1,n2)
vo_dz(m,k)=UV_ib(2,k-1,n2)
T_dz(m,k)=Tclim_ib(k-1,n2)
S_dz(m,k)=Sclim_ib(k-1,n2)
else
uo_dz(m,k)=0.5*(UV_ib(1,k-1,n2)+UV_ib(1,k,n2))
vo_dz(m,k)=0.5*(UV_ib(2,k-1,n2)+UV_ib(2,k,n2))
T_dz(m,k)=0.5*(Tclim_ib(k-1,n2)+ Tclim_ib(k,n2))
S_dz(m,k)=0.5*(Sclim_ib(k-1,n2)+ Sclim_ib(k,n2))
end if
if (k.eq.nlevels_nod2D(n2)) then ! LA 2023-08-31
uo_keel(m)=UV_ib(1,k-1,n2)
vo_keel(m)=UV_ib(2,k-1,n2)
T_keel(m)=Tclim_ib(k-1,n2)
S_keel(m)=Sclim_ib(k-1,n2)
else
uo_keel(m)=UV_ib(1,k,n2)
vo_keel(m)=UV_ib(2,k,n2)
T_keel(m)=Tclim_ib(k,n2)
S_keel(m)=Sclim_ib(k,n2)
end if
end if
end if !cavity
end do innerloop
end do nodeloop !loop over all nodes of iceberg element
do m=1,3
if (arr_ib_n_lvls(m)==0) then
cycle
else
ib_n_lvls_old = ib_n_lvls
ib_n_lvls = arr_ib_n_lvls(m)
if ((ib_n_lvls_old.ne.0) .and. ib_n_lvls_old<ib_n_lvls) ib_n_lvls=ib_n_lvls_old
end if
end do
contains
real function interpol1D(x0,f0,x1,f1,x)
implicit none
real, intent(IN) :: x0,f0,x1,f1,x
real :: frac
frac = (f1 - f0)/(x1 - x0)
interpol1D = f0 + frac * (x - x0)
end function interpol1D
end subroutine iceberg_levelwise_andkeel
!***************************************************************************************************************************
!***************************************************************************************************************************
subroutine iceberg_average_andkeel(mesh, partit, dynamics, uo_dz,vo_dz, uo_keel,vo_keel, T_dz,S_dz, T_keel,S_keel, depth_ib,iceberg_elem, ib)
USE MOD_MESH
use o_param
use MOD_PARTIT
use MOD_DYN
use o_arrays, only: Tclim_ib, Sclim_ib !, UV_ib, Z_3d_n_ib
use g_clock
use g_forcing_arrays
use g_rotate_grid
implicit none
REAL, DIMENSION(3), INTENT(OUT) :: uo_dz
REAL, DIMENSION(3), INTENT(OUT) :: vo_dz
REAL, DIMENSION(3), INTENT(OUT) :: uo_keel
REAL, DIMENSION(3), INTENT(OUT) :: vo_keel
REAL, DIMENSION(3), INTENT(OUT) :: T_dz
REAL, DIMENSION(3), INTENT(OUT) :: S_dz
REAL, DIMENSION(3), INTENT(OUT) :: T_keel
REAL, DIMENSION(3), INTENT(OUT) :: S_keel
REAl, INTENT(IN) :: depth_ib
INTEGER, INTENT(IN) :: iceberg_elem, ib
REAL, dimension(:,:,:), pointer :: UV_ib
real :: lev_up, lev_low
integer :: m, k, n2, n_up, n_low, cavity_count
! depth over which is integrated (layer and sum)
real :: dz, ufkeel1, ufkeel2, Temkeel, Salkeel
type(t_mesh), intent(in) , target :: mesh
type(t_dyn), intent(in) , target :: dynamics
type(t_partit), intent(inout), target :: partit
#include "associate_part_def.h"
#include "associate_mesh_def.h"
#include "associate_part_ass.h"
#include "associate_mesh_ass.h"
UV_IB => dynamics%uv_ib(:,:,:)
cavity_count=0
!LOOP: over all nodes of the iceberg element
nodeloop: do m=1, 3
!for each 2D node of the iceberg element..
n2=elem2D_nodes(m,iceberg_elem)
uo_dz(m)=0.0
vo_dz(m)=0.0
uo_keel(m)=0.0
vo_keel(m)=0.0
T_dz(m)=0.0
S_dz(m)=0.0
T_keel(m)=0.0
S_keel(m)=0.0
! LOOP: consider all neighboring pairs (n_up,n_low) of 3D nodes
! below n2..
!innerloop: do k=1, nl+1
innerloop: do k=1, nlevels_nod2D(n2)
lev_up = mesh%zbar_3d_n(k, n2)
if( k==nlevels_nod2D(n2) ) then
lev_low = mesh%zbar_n_bot(n2)
else
lev_low = mesh%zbar_3d_n(k+1, n2)
end if
!if( k==1 ) then
! lev_up = 0.0
!else
! lev_up = mesh%Z_3d_n_ib(k-1, n2)
! !lev_up = mesh%Z_3d_n_ib(k-1, n2)
!end if
!if( k==nlevels_nod2D(n2) ) then
! lev_low = mesh%zbar_n_bot(n2)
!else
! lev_low = mesh%Z_3d_n_ib(k, n2)
!end if
dz = abs( lev_low - lev_up )
!if( abs(lev_up)>=abs(depth_ib) ) then
! ! ...icb bottom above lev_up --> no further integration
!end if
!if( (abs(coord_nod3D(3, n_low))>abs(depth_ib)) .AND. (abs(coord_nod3D(3, n_up))>abs(depth_ib)) ) then
! write(*,*) 'INFO, k:',k,'z_up:',coord_nod3D(3, n_up),'z_lo:',coord_nod3D(3, n_low),'depth:',depth_ib,'cavity:',(mesh%cavity_flag_n(elem2D_nodes(m,iceberg_elem))==1)
!end if
if (use_cavity .AND. mesh%cavity_depth(elem2D_nodes(m,iceberg_elem)) /= 0.0 .AND. abs(depth_ib) < abs(lev_up)) then
! LA: Never go here for k=1, because abs(depth_ib)>=0.0 for all icebergs
uo_dz(m)=UV_ib(1,k-1,n2)*abs(depth_ib)
vo_dz(m)=UV_ib(2,k-1,n2)*abs(depth_ib)
uo_keel(m)=UV_ib(1,k-1,n2)
vo_keel(m)=UV_ib(2,k-1,n2)
T_dz(m)=Tclim_ib(k-1,n2)*abs(depth_ib)
S_dz(m)=Sclim_ib(k-1,n2)*abs(depth_ib)
T_keel(m)=Tclim_ib(k-1,n2)
S_keel(m)=Sclim_ib(k-1,n2) ! check those choices with RT: OK
exit innerloop
! if the lowest z coord is below the iceberg draft, exit
!else if( abs(coord_nod3D(3, n_low))>=abs(depth_ib) .AND. abs(coord_nod3D(3, n_up))<=abs(depth_ib) ) then
!****************************************************************
! LA 23.11.21 case if depth_ib<lev_up
else !# comp cav flag
if( abs(lev_low)>=abs(depth_ib) ) then !.AND. (abs(lev_up)<=abs(depth_ib)) ) then
if( abs(lev_up)<abs(depth_ib) ) then
dz = abs ( lev_up - depth_ib )
else
! LA: Never go here, when starting with k=1
dz = abs(depth_ib)
end if
!****************************************************************
if( k==1 ) then
ufkeel1 = UV_ib(1,k,n2)
ufkeel2 = UV_ib(2,k,n2)
Temkeel = Tclim_ib(k,n2)
Salkeel = Sclim_ib(k,n2)
uo_dz(m)=ufkeel1*dz
vo_dz(m)=ufkeel2*dz
T_dz(m)=Temkeel*dz
S_dz(m)=Salkeel*dz
else if( k.eq.nlevels_nod2D(n2) ) then
ufkeel1 = UV_ib(1,k-1,n2)
ufkeel2 = UV_ib(2,k-1,n2)
Temkeel = Tclim_ib(k-1,n2)
Salkeel = Sclim_ib(k-1,n2)