-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathisopyc.F
1617 lines (1591 loc) · 53.9 KB
/
isopyc.F
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
#if defined isopycmix || defined gent_mcwilliams
subroutine isopi (error, cifdef, ifdmax, nifdef, am, ah)
c
c=======================================================================
c
c Initialization for isopycnal mixing scheme
c
c -Disopycmix gives Redi/Cox version using the mixing tensor from
c Gent_Mcwilliams but no explicit use of advective velocity to
c parameterize effect of eddies on isopycnals
c
c -Disopycmix -Dgent_mcwilliams gives adds the advective velocity
c parameterization of Gent_McWilliams
c
c input:
c error = logical to signal problems
c cifdef = array of character strings for listing enabled "ifdefs"
c ifdmax = size of "cifdef"
c nifdef = current number of enabled "ifdefs"
c slmxr = 1/(max slope of isopycnals)
c ahisop = isopycnal tracer diffusivity(cm**2/sec)
c athkdf = isopycnal thickness diffusivity (cm**2/sec)
# ifdef isopycmixspatialvar
c dciso1 = isopycnal tracer diffusivity coeffs modified based
c on the slopes of the isopycnal surfaces on the east face
c of "T" cells.
c dciso2 = isopycnal tracer diffusivity coeffs modified based
c on the slopes of the isopycnal surfaces on the north face
c of "T" cells.
c dslope = half length of the interval in which "ahisop" changes
c with a steep slope from about 0.9*"ahisop" to about
c 0.1*"ahisop"
c slopec = slope at which "ahisop" is equal to half of its
c original value
# endif
c dptlim = depth limits for the reference pressure levels (in cm).
c
c output:
c The above input is potentially reset via namelist
c
c
c authors:
c R.C.Pacanowski (implemented Gokhan's version for MOM 2)
c rcp@gfdl.gov
c Gokhan Danabasoglu (pre MOM 2 version including
c Gent/Mcwilliams transport velocities)
c gokhan@isis.cgd.ucar.edu
c K. Dixon (implemented Cox version in MOM 1)
c kd@gfdl.gov
c M. Cox (original isopycnal code in COX model)
c=======================================================================
c
logical hmixset, error
character*(*) cifdef(ifdmax)
c
# include "param.h"
# include "accel.h"
# include "coord.h"
# include "grdvar.h"
# include "iounit.h"
# include "isopyc.h"
# include "levind.h"
# include "mw.h"
# include "switch.h"
# include "vmixc.h"
c
namelist /isopyc/ slmxr, ahisop, athkdf, dptlim, dslope, slopec
c
write (stdout,'(/,20x,a,/,20x,a,/)')
& 'I S O P Y C M I X I N I T I A L I Z A T I O N'
# ifdef gent_mcwilliams
&,'Gent-McWilliams version including thickness diffusion'
# else
&,'Redi/Cox version. (Gent-McWilliams without thickness diffusion)'
# endif
c
c-----------------------------------------------------------------------
c USER INPUT
c initialize variables (all mixing units are cm**2/sec.)
c-----------------------------------------------------------------------
c
slmxr = 100.0
ahisop = 2.e7
c
c define the isopycnal thickness diffusion coefficient
c
athkdf = 2.0e7
c
c reference pressure level intervals are defined (see "isopyc.h").
c "dptlim" must have "nrpl+1" elements (see "param.h"). also,
c "dptlim" are in cm.
c
c
c REMARK: the first and the last elements of "dptlim" must be the
c depth at the top (0cm) and the maximum bottom depth,
c respectively. Also, the elements of "dptlim" must be in
c increasing order.
c
dptlim(1) = 0.0e2
dptlim(2) = 1000.0e2
dptlim(3) = 2000.0e2
dptlim(4) = 3000.0e2
dptlim(5) = 4000.0e2
dptlim(6) = zw(km)
c
# ifdef isopycmixspatialvar
dslope = 0.001
slopec = 0.004
# endif
c
c-----------------------------------------------------------------------
c provide for namelist over-ride of above settings + documentation
c-----------------------------------------------------------------------
c
call getunit (io, 'namelist', 'fsr')
read (io,isopyc,end=100)
100 continue
write (stdout,isopyc)
call relunit (io)
call getunit (iodoc, 'document.dta', 'f s a')
write (iodoc, isopyc)
call relunit (iodoc)
c
c-----------------------------------------------------------------------
c add character string to "ifdef option list" indicating that this
c option is enabled
c-----------------------------------------------------------------------
c
nifdef = nifdef + 1
cifdef(nifdef) = 'isopycmix '
c
c-----------------------------------------------------------------------
c check for problems
c-----------------------------------------------------------------------
c
if (nrpl+1 .ne. 6) then
error = .true.
write (stdout,*) '=>Error: "dptlim" is set for nrpl=5 in isopi'
&, ' but nrpl=',nrpl,'. reset dptlim in isop0i or change nrpl'
endif
# if defined gent_mcwilliams && !defined isopycmix
write (stdout,*)
& '==> Error: "isopycmix" must be enabled with "gent_mcwilliams"'
error = .true.
# endif
c
c-----------------------------------------------------------------------
c print out tracer diffusion coefficients in isopycnal case
c-----------------------------------------------------------------------
c
# ifdef held_larichev
write(stdout,9102) athkdf, ah, slmxr, am
call getunit (iodoc, 'document.dta', 'f s a')
write (iodoc,'(a,e14.7)') 'am=', am, 'ah=',ah
&, 'athkdf=',athkdf, 'slmxr=',slmxr
call relunit (iodoc)
c
9102 format(/' ahisop is from held_larichev along isopyncal surfaces '
&, '(cm**2/sec) '/' athkdf = ',e12.6,' isopycnal thickness '
&,'diffusion (cm**2/sec) '/' ah = ',e12.6,' cm**2/sec for '
&,'background horizontal tracer diffusion'/' slmxr = ',e12.6
&,', to limit the slope used in computing mixing tensors',/
&,' am = ',e12.6,' cm**2/sec for horizontal mixing of momentum'/)
# else
write(stdout,9102) ahisop, athkdf, ah, slmxr, am
call getunit (iodoc, 'document.dta', 'f s a')
write (iodoc,'(a,e14.7)') 'am=', am, 'ah=',ah
&, 'ahisop=',ahisop, 'athkdf=',athkdf, 'slmxr=',slmxr
call relunit (iodoc)
c
9102 format(/' ahisop = ',e12.6,' along isopyncal tracer mixing '
&, '(cm**2/sec) '/' athkdf = ',e12.6,' isopycnal thickness '
&,'diffusion (cm**2/sec) '/' ah = ',e12.6,' cm**2/sec for '
&,'background horizontal tracer diffusion'/' slmxr = ',e12.6
&,', to limit the slope used in computing mixing tensors',/
&,' am = ',e12.6,' cm**2/sec for horizontal mixing of momentum'/)
# endif
c
c-----------------------------------------------------------------------
c store the square root of the tracer timestep acceleration values
c into variable "dtxsqr" for use in isopycnal mixing
c-----------------------------------------------------------------------
c
do k=1,km
dtxsqr(k) = sqrt(dtxcel(k))
enddo
c
c-----------------------------------------------------------------------
c determine the isopycnal reference pressure levels for the "t"
c grid point levels, using the depths at the "t" grid points as the
c reference depth (pressure)
c-----------------------------------------------------------------------
c
do k=1,km
do m=2,nrpl+1
if (zt(k).gt.dptlim(m-1) .and. zt(k).le.dptlim(m)) then
kisrpl(k) = m-1
go to 101
endif
enddo
101 continue
if (kisrpl(k) .lt. 1 .or. kisrpl(k) .gt. nrpl) then
write (stdout,9100) kisrpl(k), k
error = .true.
endif
enddo
9100 format (/,' =>Error: kisrpl is ',i3,' at k ',1x,i3,' in isopyc.F')
c
c-----------------------------------------------------------------------
c the indices used in isopycnal mixing indicating the location of
c the reference pressure levels in the 20-level table of polynomial
c expansion variables are computed
c
c REMARK: because the polynomial expansion coefficients are
c functions of the reference potential temperature and
c salinity profiles, at the reference pressure level
c the corresponding potential temperature and salinity
c values will be used.
c-----------------------------------------------------------------------
c
do m=1,nrpl
krplin(m) = 0
enddo
c
do m=2,nrpl+1
dptmid = p5*(dptlim(m-1)+dptlim(m))
if (dptmid .le. zt(1)) then
krplin(m-1) = 1
elseif (dptmid .gt. zt(km)) then
krplin(m-1) = km
elseif (dptmid.gt.zt(1) .and. dptmid.le.zt(km)) then
do k=2,km
if (zt(k) .ge. dptmid) then
t1 = zt(k)-dptmid
t2 = dptmid-zt(k-1)
if (t1 .gt. t2) then
krplin(m-1) = k-1
else
krplin(m-1) = k
endif
go to 102
endif
enddo
102 continue
endif
if (krplin(m-1) .lt. 1 .or. krplin(m-1) .gt. km) then
write (stderr,9110) krplin(m-1),m-1
error = .true.
endif
enddo
c
write (stdout,96) (kisrpl(k),k=1,km)
write (stdout,97) (krplin(m),m=1,nrpl)
c
c-----------------------------------------------------------------------
c the isopycnal diffusion coefficient may be a function of depth. in
c the default configuration, the isopycnal diffusion coefficient is
c a constant: "fzisop", which multiplies "ahisop", is set to unity.
c if "ahisop" varies in the vertical, "fzisop" should contain this
c variation. the value of "ahisop" should be adjusted accordingly.
c-----------------------------------------------------------------------
c
do k=1,km
fzisop(k) = c1
enddo
c
write (stdout,98)
write (stdout,99) (fzisop(k),k=1,km)
c
96 format (/,' isopycnal reference pressure levels (kisrpl) = ',
& 20(1x,i4))
97 format (/,' reference pressure level indices (krplin) = ',
& 20(1x,i4))
98 format (/,' vertical variation of "ahisop" (fzisop) = ')
99 format (5(1x,e12.6))
9110 format (/,' =>Error: krplin is ',i3,' at m ',1x,i3)
c
c-----------------------------------------------------------------------
c initialize arrays
c-----------------------------------------------------------------------
c
do j=jsmw,jemw
do i=1,imt
do k=1,km
K1(i,k,j,3) = c0
K3(i,k,j,1) = c0
K3(i,k,j,2) = c0
K3(i,k,j,3) = c0
# ifdef gent_mcwilliams
adv_vetiso(i,k,j) = c0
# ifdef isopycmixspatialvar
dciso1(i,k,j) = c0
# endif
# endif
enddo
enddo
enddo
c
do j=1,jemw
do i=1,imt
do k=1,km
K2(i,k,j,3) = c0
# ifdef gent_mcwilliams
adv_vntiso(i,k,j) = c0
# ifdef isopycmixspatialvar
dciso2(i,k,j) = c0
# endif
# endif
enddo
enddo
enddo
c
do m=1,nrpl
do j=1,jmw
do i=1,imt
do k=1,km
rhoi(i,k,j,m) = c0
enddo
enddo
enddo
enddo
c
do m=1,3
do j=1,jemw
do i=1,imt
do k=1,km+1
e(i,k,j,m) = c0
enddo
enddo
enddo
enddo
c
# ifdef gent_mcwilliams
do j=jsmw,jemw
do k=0,km
do i=1,imt
adv_vbtiso(i,k,j) = c0
adv_fbiso(i,k,j) = c0
enddo
enddo
enddo
# endif
return
end
subroutine isopyc (joff, js, je, is, ie)
c
c=======================================================================
c
c Compute the isopycnal mixing tensor components and the
c isopycnal advection velocities which parameterize the effect
c of eddies on the isopycnals.
c
c
c Mixing tensor "K" is ...
c
c | 1.0 K1(,,,2) K1(,,,3) |
c | |
c K = | K2(,,,1) 1.0 K2(,,,3) |
c | |
c | K3(,,,1) K3(,,,2) K3(,,,3) |
c
c where K1(,,,2) and K2(,,,1) are set to 0.0 (neglected)
c
c
c input:
c joff = offset relating row "j" in the MW to latitude "jrow"
c js = starting row within the MW for calculations
c je = ending row within the MW for calculations
c is = starting index longitude within the MW
c ie = ending index longitude within the MW
c
c output:
c rhoi = density at tau-1 referenced to pressure levels
c K1 = tensor components (1,2), (1,3) centered on east face
c of "T" cells
c K2 = tensor components (2,1), (2,3) centered on north face
c of "T" cells
c K3 = tensor components (3,1), (3,2), (3,3) centered on
c bottom face of "T" cells
# ifdef gent_mcwilliams
c adv_vetiso = isopycnal advective vel on east face of "T" cell
c adv_vntiso = isopycnal advective vel on north face of "T" cell
c (Note: this includes the cosine factor as in "adv_vnt")
c adv_vbtiso = isopycnal advective vel on bottom face of "T" cell
# endif
c
c=======================================================================
c
# include "param.h"
# include "accel.h"
# include "grdvar.h"
# include "iounit.h"
# include "isopyc.h"
# include "levind.h"
# include "mw.h"
# include "switch.h"
# include "state.h"
# include "dens.h"
c
# ifdef timing
call tic ('isopycnal ', 'isopyc')
# endif
c
c-----------------------------------------------------------------------
c set local constants
c-----------------------------------------------------------------------
c
istrt = max(2,is)
iend = min(imt-1,ie)
c
c-----------------------------------------------------------------------
c after the first MW, move variables from the top two rows to the
c bottom two rows of the MW to eliminate redundant calculation
c-----------------------------------------------------------------------
c
if (joff .ne. 0) then
call moveiso (istrt-1, iend+1)
endif
c
c-----------------------------------------------------------------------
c compute normalized densities for each isopycnal reference pressure
c level using a 3rd order polynomial fit to the equation of state.
c for each isopycnal reference pressure level, the same reference
c potential temperature, reference salinity and expansion coeff
c values are used at all of the vertical levels.
c
c Note: this density is used for the mixing tensor in both the
c Redi/Cox and Gent/McWilliams options
c
c first MW : calculate for rows 1..jmw
c 1 < MW < last: calculate for rows jsmw+1..jmw
c last MW : calculate for rows jsmw+1..<=jmw
c-----------------------------------------------------------------------
c
do m=1,nrpl
tref = to(krplin(m))
sref = so(krplin(m))
do j=js, je
do k=1,km
do i=1,imt
tprime = t(i,k,j,1,taum1) - tref
sprime = t(i,k,j,2,taum1) - sref
rhoi(i,k,j,m) = dens (tprime, sprime, krplin(m))
enddo
enddo
enddo
enddo
#ifdef trace_indices
write (stdout,'(2x,5(a,i4))')
& "=> In isopyc: rhoi from js=",js," je=",je," joff=",joff
&," jrows=",js+joff," to ",je+joff
#endif
c
c-----------------------------------------------------------------------
c evaluate K2(,,3) centered on the northern face of "T" cells
c first MW : calculate for rows 1..jemw
c 1 < MW < last: calculate for rows jsmw..jemw
c last MW : calculate for rows jsmw..<=jemw
c-----------------------------------------------------------------------
c
call k2_3 (joff, max(js-1,1), je-1, is, ie)
c
c-----------------------------------------------------------------------
c evaluate K1(,,3) centered on eastern face of "T" cells
c 1 <= MW < last: calculate for rows jsmw..jemw
c last MW : calculate for rows jsmw..<=jemw
c-----------------------------------------------------------------------
c
call k1_3 (joff, max(js-1,jsmw), je-1, is, ie)
c
c-----------------------------------------------------------------------
c evaluate K3(,,1..3) centered on bottom face of "T" cells
c 1 <= MW < last: calculate for rows jsmw..jemw
c last MW : calculate for rows jsmw..<=jemw
c-----------------------------------------------------------------------
c
call k3_123 (joff, max(js-1,jsmw), je-1, is, ie)
c
# ifdef gent_mcwilliams
c
c-----------------------------------------------------------------------
c compute isopycnal advective velocities for tracers
c first MW : calculate (viso) for rows 1..jemw
c (uiso,wiso) for rows 2..jemw
c 1 < MW < last: calculate (uiso,viso,wiso) for rows jsmw..jemw
c last MW : calculate (uiso,viso,wiso) for rows jsmw..<=jemw
c-----------------------------------------------------------------------
c
call isopyc_adv (joff, max(js-1,1), je-1, is, ie)
c
# endif
c
# ifdef timing
call toc ('isopycnal ', 'isopyc')
# endif
return
end
subroutine moveiso (is, ie)
c
c=======================================================================
c move the MW up (northward) by moving data from the last two rows
c into the first two rows.
c
c input:
c is = starting longitude index in the MW
c ie = ending longitude index in the MW
c=======================================================================
c
#include "param.h"
#include "hmixc.h"
#include "isopyc.h"
c
nrows = jmw - ncrows
do move=1,nrows
jfrom = jmw - (nrows - move)
jto = move
c
c-----------------------------------------------------------------------
c move quantities with rows dimensioned (1:jmw)
c-----------------------------------------------------------------------
c
do n=1,nrpl
do k=1,km
do i=is,ie
rhoi(i,k,jto,n) = rhoi(i,k,jfrom,n)
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c move quantities with rows dimensioned (1:jemw)
c-----------------------------------------------------------------------
c
if (jfrom .le. jemw) then
do k=1,km
do i=is,ie
K2(i,k,jto,3) = K2(i,k,jfrom,3)
# ifdef gent_mcwilliams
adv_vntiso(i,k,jto) = adv_vntiso(i,k,jfrom)
# endif
enddo
enddo
endif
# ifdef trace_indices
write (stdout,'(4x,2(a,i4))')
& "=> In moveiso: moving variables on row ",jfrom," to row ",jto
# endif
enddo
c
return
end
subroutine k1_3 (joff, js, je, is, ie)
c
c=======================================================================
c
c compute "K1(,,3)" at the center of the eastern face of "T" cells
c use "c1e10" to keep the exponents in range.
c
c=======================================================================
c
# include "param.h"
# include "accel.h"
# include "coord.h"
# include "grdvar.h"
# include "iounit.h"
# include "isopyc.h"
# include "levind.h"
# include "mw.h"
# include "switch.h"
c
c-----------------------------------------------------------------------
c set local constants
c-----------------------------------------------------------------------
c
c1e10 = 1.0e10
eps = 1.0e-25
c
c-----------------------------------------------------------------------
c d(rho_barx_barz)/dz centered on eastern face of "t" cells
c Note: values involving ocean surface and ocean bottom are
c estimated afterwards using a linear extrapolation
c-----------------------------------------------------------------------
c
do j=js,je
do k=2,kmm1
m = kisrpl(k)
fxd = c1e10*p5*dzt2r(k)
do i=1,imtm1
e(i,k,j,3) = fxd*(rhoi(i ,k-1,j,m) - rhoi(i ,k+1,j,m)
& +rhoi(i+1,k-1,j,m) - rhoi(i+1,k+1,j,m))
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c linearly extrapolate densities to ocean surface for calculation
c of d(rho_barx_barz)/dz involving level 1.
c
c REMARK: requires min(kmt(i,jrow)) = 2 cells in ocean.
c-----------------------------------------------------------------------
c
k = 1
fxd = c1e10*dztr(k)
fxe = dzw(k-1)+dzw(k)
m = kisrpl(k)
do j=js,je
do i=1,imtm1
fxa = p5*(rhoi(i,k+1,j,m) + rhoi(i+1,k+1,j,m))
fxb = p5*(rhoi(i,k,j,m) + rhoi(i+1,k,j,m))
fxc = dzwr(k)*(fxb*fxe - fxa*dzw(k-1))
e(i,k,j,3) = fxd*(fxc - p5*(fxa+fxb))
enddo
enddo
c
c-----------------------------------------------------------------------
c linearly extrapolate densities to ocean bottom for calculation
c of d(rho_barx_barz)/dz involving bottom level.
c-----------------------------------------------------------------------
c
do j=js,je
do i=1,imtm1
e(i,km,j,3) = c0
enddo
enddo
c
do j=js,je
jrow = j + joff
do i=1,imtm1
k = min(kmt(i,jrow),kmt(i+1,jrow))
if (k .ne. 0) then
fxe = dzw(k-1)+dzw(k)
m = kisrpl(k)
fxa = p5*(rhoi(i,k-1,j,m) + rhoi(i+1,k-1,j,m))
fxb = p5*(rhoi(i,k ,j,m) + rhoi(i+1,k ,j,m))
fxc = dzwr(k-1)*(fxb*fxe - fxa*dzw(k))
e(i,k,j,3) = dztr(k)*c1e10*(p5*(fxa+fxb) - fxc)
endif
enddo
enddo
c
c-----------------------------------------------------------------------
c "e(,,,1)" = d(rho)/dx centered on east face of "T" cells
c "e(,,,2)" = d(rho_barx_bary)/dy centered on east face of "T" cells
c-----------------------------------------------------------------------
c
do j=js,je
jrow = j + joff
do k=1,km
m = kisrpl(k)
do i=1,imtm1
e(i,k,j,1) = tmask(i,k,j)*tmask(i+1,k,j)*cstr(jrow)*dxur(i)
& *c1e10*(rhoi(i+1,k,j,m) - rhoi(i,k,j,m))
e(i,k,j,2) = dyt4r(jrow)*c1e10*(
& rhoi(i ,k,j+1,m) - rhoi(i ,k,j-1,m)
& + rhoi(i+1,k,j+1,m) - rhoi(i+1,k,j-1,m))
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c if any one of the 4 neighboring corner grid points is a land point,
c set "e(i,k,j,2)" to zero. note that "e(i,k,j,2)" will be used
c only in the slope check.
c-----------------------------------------------------------------------
c
do j=js,je
do k=1,km
do i=1,imtm1
olmask = tmask(i,k,j-1)*tmask(i,k,j+1)*tmask(i+1,k,j-1)
& *tmask(i+1,k,j+1)
if (olmask .eq. c0) e(i,k,j,2) = c0
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c impose zonal boundary conditions at "i"=1 and "imt"
c-----------------------------------------------------------------------
c
do j=js,je
call setbcx (e(1,1,j,1), imt, km)
call setbcx (e(1,1,j,2), imt, km)
call setbcx (e(1,1,j,3), imt, km)
enddo
c
c-----------------------------------------------------------------------
c compute "K1", using "slmxr" to limit vertical slope of isopycnal
c to guard against numerical instabilities.
c-----------------------------------------------------------------------
c
do j=js,je
do k=1,km
# ifdef isopycmixspatialvar
do i=1,imt
fxa = c0
fxb = sign(1.,e(i,k,j,3))/(abs(e(i,k,j,3))+eps)
slope = fxb*sqrt(e(i,k,j,1)**2+e(i,k,j,2)**2)
if (slope .le. c0 .and. slope .ge. (-c1/slmxr)) then
fxa = p5*(c1+tanh((slope+slopec)/dslope))
endif
dciso1(i,k,j) = fxa
K1(i,k,j,3) = -fxb*e(i,k,j,1)*dciso1(i,k,j)
enddo
# else
do i=1,imt
chkslp = -sqrt(e(i,k,j,1)**2+e(i,k,j,2)**2)*slmxr*dtxsqr(k)
if (e(i,k,j,3) .gt. chkslp) e(i,k,j,3) = chkslp
enddo
do i=1,imt
K1(i,k,j,3) = (-e(i,k,j,1)*e(i,k,j,3)*fzisop(k))
& /(e(i,k,j,3)**2+eps)
enddo
# endif
enddo
enddo
#ifdef trace_indices
write (stdout,'(2x,5(a,i4))')
& "=> In K1_3: js=",js," je=",je," joff=",joff
&," jrows=",js+joff," to ",je+joff
#endif
return
end
subroutine k2_3 (joff, js, je, is, ie)
c
c=======================================================================
c compute "K2(,,3)" at the center of the northern face of "T" cells
c use "c1e10" to keep the exponents in range.
c=======================================================================
c
# include "param.h"
# include "accel.h"
# include "coord.h"
# include "grdvar.h"
# include "iounit.h"
# include "isopyc.h"
# include "levind.h"
# include "mw.h"
# include "switch.h"
c
c-----------------------------------------------------------------------
c set local constants
c-----------------------------------------------------------------------
c
c1e10 = 1.0e10
eps = 1.0e-25
c
c-----------------------------------------------------------------------
c d(rho_bary_barz)/dz centered on northern face of "T" cells
c Note: values involving ocean surface and ocean bottom are
c estimated afterwards using a linear extrapolation
c-----------------------------------------------------------------------
c
do j=js,je
do k=2,kmm1
m = kisrpl(k)
fxd = c1e10*p5*dzt2r(k)
do i=2,imtm1
e(i,k,j,3) = fxd*(rhoi(i,k-1,j ,m) - rhoi(i,k+1,j ,m)
& +rhoi(i,k-1,j+1,m) - rhoi(i,k+1,j+1,m))
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c linearly extrapolate densities to ocean surface for calculation
c of d(rho_bary_barz)/dz involving level 1.
c
c REMARK: requires min(kmt(i,jrow)) = 2 cells in ocean.
c-----------------------------------------------------------------------
c
k = 1
fxd = c1e10*dztr(k)
fxe = dzw(k-1)+dzw(k)
m = kisrpl(k)
do j=js,je
do i=2,imtm1
fxa = p5*(rhoi(i,k+1,j,m) + rhoi(i,k+1,j+1,m))
fxb = p5*(rhoi(i,k ,j,m) + rhoi(i,k ,j+1,m))
fxc = dzwr(k)*(fxb*fxe - fxa*dzw(k-1))
e(i,k,j,3) = fxd*(fxc - p5*(fxa+fxb))
enddo
enddo
c
c-----------------------------------------------------------------------
c linearly extrapolate densities to ocean bottom for calculation
c of d(rho_bary_barz)/dz involving bottom level.
c-----------------------------------------------------------------------
c
do j=js,je
do i=2,imtm1
e(i,km,j,3) = c0
enddo
enddo
c
do j=js,je
jrow = j + joff
do i=2,imtm1
k = min(kmt(i,jrow),kmt(i,jrow+1))
if (k .ne. 0) then
fxe = dzw(k-1)+dzw(k)
m = kisrpl(k)
fxa = p5*(rhoi(i,k-1,j,m) + rhoi(i,k-1,j+1,m))
fxb = p5*(rhoi(i,k ,j,m) + rhoi(i,k ,j+1,m))
fxc = dzwr(k-1)*(fxb*fxe - fxa*dzw(k))
e(i,k,j,3) = dztr(k)*c1e10*(p5*(fxa+fxb) - fxc)
endif
enddo
enddo
c
c-----------------------------------------------------------------------
c "e(,,,1)" = d(rho_barx_bary)/dx centered on north face of "T" cells
c "e(,,,2)" = d(rho)/dy on north face of "T" cells
c-----------------------------------------------------------------------
c
do j=js,je
jrow = j + joff
do k=1,km
m = kisrpl(k)
do i=2,imtm1
e(i,k,j,1) = csur(jrow)*dxt4r(i)*c1e10*(
& rhoi(i+1,k,j+1,m) - rhoi(i-1,k,j+1,m)
& + rhoi(i+1,k,j ,m) - rhoi(i-1,k,j ,m))
e(i,k,j,2) = tmask(i,k,j)*tmask(i,k,j+1)*dyur(jrow)*c1e10
& *(rhoi(i,k,j+1,m) - rhoi(i,k,j,m))
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c if any one of the 4 neighboring corner grid points is a land point,
c set "e(i,k,j,1)" to zero. note that "e(i,k,j,1)" will be used
c only in the slope check.
c-----------------------------------------------------------------------
c
do j=js,je
do k=1,km
do i=2,imtm1
olmask = tmask(i-1,k,j+1)*tmask(i+1,k,j+1)*tmask(i-1,k,j)
& *tmask(i+1,k,j)
if (olmask .eq. c0) e(i,k,j,1) = c0
enddo
enddo
enddo
c
c-----------------------------------------------------------------------
c impose zonal boundary conditions at "i"=1 and "imt"
c-----------------------------------------------------------------------
c
do j=js,je
call setbcx (e(1,1,j,1), imt, km)
call setbcx (e(1,1,j,2), imt, km)
call setbcx (e(1,1,j,3), imt, km)
enddo
c
c-----------------------------------------------------------------------
c compute "K2", using "slmxr" to limit vertical slope of isopycnal
c to guard against numerical instabilities.
c-----------------------------------------------------------------------
c
do j=js,je
do k=1,km
# ifdef isopycmixspatialvar
do i=1,imt
fxa = c0
fxb = sign(1.,e(i,k,j,3))/(abs(e(i,k,j,3))+eps)
slope = fxb*sqrt(e(i,k,j,1)**2+e(i,k,j,2)**2)
if (slope .le. c0 .and. slope .ge. (-c1/slmxr)) then
fxa = p5*(c1+tanh((slope+slopec)/dslope))
endif
dciso2(i,k,j) = fxa
K2(i,k,j,3) = -fxb*e(i,k,j,2)*dciso2(i,k,j)
enddo
# else
do i=1,imt
chkslp = -sqrt(e(i,k,j,1)**2+e(i,k,j,2)**2)*slmxr*dtxsqr(k)
if (e(i,k,j,3) .gt. chkslp) e(i,k,j,3) = chkslp
enddo
do i=1,imt
K2(i,k,j,3) = (-e(i,k,j,2)*e(i,k,j,3)*fzisop(k))
& /(e(i,k,j,3)**2+eps)
enddo
# endif
enddo
enddo
#ifdef trace_indices
write (stdout,'(2x,5(a,i4))')
& "=> In K2_3: js=",js," je=",je," joff=",joff
&," jrows=",js+joff," to ",je+joff
#endif
return
end
subroutine k3_123 (joff, js, je, is, ie)
c
c=======================================================================
c compute K2(,,,1:3) at the center of the bottom face of "T" cells
c use "c1e10" to keep the exponents in range.
c=======================================================================
c
# include "param.h"
# include "accel.h"
# include "grdvar.h"
# include "iounit.h"
# include "isopyc.h"
# include "levind.h"
# include "mw.h"
# include "switch.h"
# include "vmixc.h"
c
c-----------------------------------------------------------------------
c set local constants
c-----------------------------------------------------------------------
c
c1e10 = 1.0e10
eps = 1.0e-25
c
do j=js,je
jrow = j + joff
do k=2,km
m = kisrpl(k)
do i=2,imtm1
e(i,k-1,j,1) = cstr(jrow)*dxt4r(i)*c1e10
& *(tmask(i-1,k-1,j)*tmask(i,k-1,j)*(rhoi(i,k-1,j,m)
& -rhoi(i-1,k-1,j,m))
& +tmask(i,k-1,j)*tmask(i+1,k-1,j)*(rhoi(i+1,k-1,j,m)
& -rhoi(i,k-1,j,m))
& +tmask(i-1,k,j)*tmask(i,k,j)*(rhoi(i,k,j,m)
& -rhoi(i-1,k,j,m))
& +tmask(i,k,j)*tmask(i+1,k,j)*(rhoi(i+1,k,j,m)
& -rhoi(i,k,j,m)))
e(i,k-1,j,2) = dyt4r(jrow)*c1e10
& *(tmask(i,k-1,j-1)*tmask(i,k-1,j)*(rhoi(i,k-1,j,m)
& -rhoi(i,k-1,j-1,m))
& +tmask(i,k-1,j)*tmask(i,k-1,j+1)*(rhoi(i,k-1,j+1,m)
& -rhoi(i,k-1,j,m))
& +tmask(i,k,j-1)*tmask(i,k,j)*(rhoi(i,k,j,m)
& -rhoi(i,k,j-1,m))
& +tmask(i,k,j)*tmask(i,k,j+1)*(rhoi(i,k,j+1,m)
& -rhoi(i,k,j,m)))
e(i,k-1,j,3) = dzwr(k-1)*tmask(i,k-1,j)*tmask(i,k,j)*c1e10
& *(rhoi(i,k-1,j,m) - rhoi(i,k,j,m))
enddo
enddo
k = km
do i=2,imtm1
e(i,k,j,1) = c0
e(i,k,j,2) = c0
e(i,k,j,3) = c0
enddo
enddo
c
c-----------------------------------------------------------------------
c compute "K3", using "slmxr" to limit vertical slope of isopycnal
c to guard against numerical instabilities.
c-----------------------------------------------------------------------
c
do j=js,je
do k=1,km
# ifdef isopycmixspatialvar
do i=2,imt-1
fxa = c0
fxb = sign(1.,e(i,k,j,3))/(abs(e(i,k,j,3))+eps)
slope = fxb*sqrt(e(i,k,j,1)**2+e(i,k,j,2)**2)
if (slope .le. c0 .and. slope .ge. (-c1/slmxr)) then
fxa = p5*(c1+tanh((slope+slopec)/dslope))
endif
fxc = fxb*fxa
K3(i,k,j,1) = -fxc*e(i,k,j,1)
K3(i,k,j,2) = -fxc*e(i,k,j,2)
K3(i,k,j,3) = fxb*fxb*(e(i,k,j,1)**2+e(i,k,j,2)**2)
& *fxa
enddo
# else
do i=2,imt-1
chkslp = -sqrt(e(i,k,j,1)**2+e(i,k,j,2)**2)*slmxr*dtxsqr(k)
if (e(i,k,j,3) .gt. chkslp) e(i,k,j,3) = chkslp
enddo
fxa = p5*(fzisop(min(km,k+1))+fzisop(k))
do i=2,imt-1
ahfctr = fxa/(e(i,k,j,3)**2+eps)
K3(i,k,j,1) = -e(i,k,j,3)*e(i,k,j,1)*ahfctr
K3(i,k,j,2) = -e(i,k,j,3)*e(i,k,j,2)*ahfctr
K3(i,k,j,3) = (e(i,k,j,1)**2+e(i,k,j,2)**2)*ahfctr