forked from HYCOM/HYCOM-src
-
Notifications
You must be signed in to change notification settings - Fork 0
/
geopar.F90
1067 lines (1066 loc) · 33.5 KB
/
geopar.F90
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
subroutine geopar
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
use mod_za ! HYCOM I/O interface
!
! --- set up model parameters related to geography
!
! --- hycom version 2.1
!
implicit none
!
real*4, allocatable, dimension(:,:) :: g_sc
!
real dx0kf,dp0kf,dpm,dpms,ds0kf,dsm,dsms
real hmina,hminb,hmaxa,hmaxb
real*8 sum_ip,sum_is,sum_isa
integer i,ios,j,k,ktr,l,nishlf
character preambl(5)*79,cline*80
!
real aspmax
parameter (aspmax=2.0) ! maximum grid aspect ratio for diffusion
! parameter (aspmax=1.0) ! ignore grid aspect ratio in diffusion
!
! --- read grid location,spacing,coriolis arrays
!
if (mnproc.eq.1) then ! .b file from 1st tile only
write (lp,'(3a)') ' reading grid file from ', &
trim(flnmgrd),'.[ab]'
open (unit=uoff+9,file=trim(flnmgrd)//'.b', &
status='old')
endif
call xcsync(flush_lp)
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
read(cline,*) i
!
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
read (cline,*) j
!
if (i.ne.itdm .or. j.ne.jtdm) then
if (mnproc.eq.1) then
write(lp,'(/ a /)') &
'error - wrong array size in grid file'
endif
call xcstop('(geopar)')
stop '(geopar)'
endif
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
if (mnproc.eq.1) then
write (lp,'(a)') trim(cline)
endif
read (cline,*) mapflg
!
call zaiopf(trim(flnmgrd)//'.a','old', 9)
!
do k= 1,15
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
i = index(cline,'=')
read (cline(i+1:),*) hminb,hmaxb
if (mnproc.eq.1) then
write (lp,'(a)') trim(cline)
endif
call xcsync(flush_lp)
!
if (k.eq.1) then
call zaiord(plon, ip,.false., hmina,hmaxa, 9)
elseif (k.eq.2) then
call zaiord(plat, ip,.false., hmina,hmaxa, 9)
do i= 1,2 !skip qlon,qlat
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ', &
uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
#if defined (USE_NUOPC_CESMBETA)
! qlon, qlat
j = index(cline,'=')
read (cline(j+1:),*) hminb,hmaxb
if (mnproc.eq.1) then
write (lp,'(a)') trim(cline)
endif
call xcsync(flush_lp)
call zaiord(util2, ip,.false., hmina,hmaxa, 9)
if (i.eq.1) then
qlon(:,:)=util2(:,:)
else
qlat(:,:)=util2(:,:)
endif
#else
! skip qlon, qlat
call zaiosk(9)
#endif
enddo
elseif (k.eq.3) then
call zaiord(ulon, ip,.false., hmina,hmaxa, 9)
elseif (k.eq.4) then
call zaiord(ulat, ip,.false., hmina,hmaxa, 9)
elseif (k.eq.5) then
call zaiord(vlon, ip,.false., hmina,hmaxa, 9)
elseif (k.eq.6) then
call zaiord(vlat, ip,.false., hmina,hmaxa, 9)
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
#if defined(ESPC_COUPLE) || defined (USE_NUOPC_CESMBETA)
! pang
i = index(cline,'=')
read (cline(i+1:),*) hminb,hmaxb
if (mnproc.eq.1) then
write (lp,'(a)') trim(cline)
endif
call xcsync(flush_lp)
call zaiord(pang, ip,.false., hmina,hmaxa, 9)
#else
! skip pang
call zaiosk(9)
#endif
elseif (k.eq.7) then
call zaiord(scpx, ip,.false., hmina,hmaxa, 9)
elseif (k.eq.8) then
call zaiord(scpy, ip,.false., hmina,hmaxa, 9)
elseif (k.eq.9) then
call zaiord(scqx, iq,.false., hmina,hmaxa, 9)
elseif (k.eq.10) then
call zaiord(scqy, iq,.false., hmina,hmaxa, 9)
elseif (k.eq.11) then
call zaiord(scux, iu,.false., hmina,hmaxa, 9)
elseif (k.eq.12) then
call zaiord(scuy, iu,.false., hmina,hmaxa, 9)
elseif (k.eq.13) then
call zaiord(scvx, iv,.false., hmina,hmaxa, 9)
elseif (k.eq.14) then
call zaiord(scvy, iv,.false., hmina,hmaxa, 9)
else
call zaiord(corio,iq,.false., hmina,hmaxa, 9)
endif
!
if (abs(hmina-hminb).gt.abs(hminb)*1.e-4 .or. &
abs(hmaxa-hmaxb).gt.abs(hmaxb)*1.e-4 ) then
if (mnproc.eq.1) then
write(lp,'(/ a / a,1p3e14.6 / a,1p3e14.6 /)') &
'error - .a and .b files not consistent:', &
'.a,.b min = ',hmina,hminb,hmina-hminb, &
'.a,.b max = ',hmaxa,hmaxb,hmaxa-hmaxb
endif
call xcstop('(geopar)')
stop '(geopar)'
endif
enddo
!
call zaiocl(9)
if (mnproc.eq.1) then ! .b file from 1st tile only
close(unit=uoff+9)
endif
!
if (itest.gt.0 .and. jtest.gt.0) then
i=itest
j=jtest
write (lp,'(/ a,2i5,a,f8.3,a,f12.9,2f10.2/)') &
' i,j=',i+i0,j+j0, &
' plat=',plat(i,j), &
' corio,scux,vy=',corio(i,j),scux(i,j),scvy(i,j)
endif
call xcsync(flush_lp)
!
! --- read basin depth array
!
if (mnproc.eq.1) then ! .b file from 1st tile only
write (lp,'(3a)') ' reading bathymetry file from ', &
trim(flnmdep),'.[ab]'
open (unit=uoff+9,file=trim(flnmdep)//'.b', &
status='old')
read ( uoff+9,'(a79)') preambl
endif
call xcsync(flush_lp)
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
i = index(cline,'=')
read (cline(i+1:),*) hminb,hmaxb
if (mnproc.eq.1) then ! .b file from 1st tile only
close(unit=uoff+9)
write (lp,'(/(1x,a))') preambl,cline
endif
!
call zaiopf(trim(flnmdep)//'.a','old', 9)
call zaiord(depths,ip,.false., hmina,hmaxa, 9)
call zaiocl(9)
!
if (abs(hmina-hminb).gt.abs(hminb)*1.e-4 .or. &
abs(hmaxa-hmaxb).gt.abs(hmaxb)*1.e-4 ) then
if (mnproc.eq.1) then
write(lp,'(/ a / a,1p3e14.6 / a,1p3e14.6 /)') &
'error - .a and .b files not consistent:', &
'.a,.b min = ',hmina,hminb,hmina-hminb, &
'.a,.b max = ',hmaxa,hmaxb,hmaxa-hmaxb
endif
call xcstop('(geopar)')
stop '(geopar)'
endif
!
!$OMP PARALLEL DO PRIVATE(j,i) &
!$OMP SCHEDULE(STATIC,jblk)
do j= 1,jj
do i= 1,ii
if (depths(i,j).gt.0.5*hugel) then
depths(i,j) = 0.0
endif
enddo
enddo
!
! --- determine do-loop limits for u,v,p,q points, and update halo for depths
call bigrid(depths, mapflg, util1,util2,util3)
!cc call prtmsk(ip,depths,util1,idm,ii,jj,0.0,1.0,
!cc & 'bottom depth (m)')
!
! now safe to apply halo to arrays.
!
vland = 1.0
call xctilr(plon, 1,1, nbdy,nbdy, halo_ps)
call xctilr(plat, 1,1, nbdy,nbdy, halo_ps)
#if defined(USE_NUOPC_CESMBETA)
call xctilr(qlon, 1,1, nbdy,nbdy, halo_ps)
call xctilr(qlat, 1,1, nbdy,nbdy, halo_ps)
#endif
#if defined(ESPC_COUPLE) || defined (USE_NUOPC_CESMBETA)
call xctilr(pang, 1,1, nbdy,nbdy, halo_ps)
#endif
if (momtyp.eq.2) then
call xctilr(scpx, 1,1, nbdy,nbdy, halo_ps)
call xctilr(scpy, 1,1, nbdy,nbdy, halo_ps)
call xctilr(scqx, 1,1, nbdy,nbdy, halo_qs)
call xctilr(scqy, 1,1, nbdy,nbdy, halo_qs)
endif !momtyp=2
call xctilr(corio, 1,1, nbdy,nbdy, halo_qs)
call xctilr(ulon, 1,1, nbdy,nbdy, halo_us)
call xctilr(ulat, 1,1, nbdy,nbdy, halo_us)
call xctilr(scux, 1,1, nbdy,nbdy, halo_us)
call xctilr(scuy, 1,1, nbdy,nbdy, halo_us)
call xctilr(vlon, 1,1, nbdy,nbdy, halo_vs)
call xctilr(vlat, 1,1, nbdy,nbdy, halo_vs)
call xctilr(scvx, 1,1, nbdy,nbdy, halo_vs)
call xctilr(scvy, 1,1, nbdy,nbdy, halo_vs)
vland = 0.0
!
! --- momtum4 needs sc[pq][xy] defined everywhere
!
if (momtyp.ne.2) then
allocate( g_sc(itdm,jtdm) )
call zaiopf(trim(flnmgrd)//'.a','old', 9)
do k= 1,9
call zaiosk(9)
enddo
call zaiordg(g_sc, 9)
call geopar_halo(scpx, g_sc, halo_ps)
call zaiordg(g_sc, 9)
call geopar_halo(scpy, g_sc, halo_ps)
call zaiordg(g_sc, 9)
call geopar_halo(scqx, g_sc, halo_qs)
call zaiordg(g_sc, 9)
call geopar_halo(scqy, g_sc, halo_qs)
call zaiocl(9)
deallocate(g_sc)
endif !momtyp/=2
!
! --- area of grid cells (length x width) at u,v,p,q points resp.
!
!*****!$OMP PARALLEL DO PRIVATE(j,i)
!*****!$OMP& SCHEDULE(STATIC,jblk)
do j=1-nbdy,jj+nbdy
do i=1-nbdy,ii+nbdy
scu2(i,j)=scux(i,j)*scuy(i,j)
scv2(i,j)=scvx(i,j)*scvy(i,j)
scp2(i,j)=scpx(i,j)*scpy(i,j)
scq2(i,j)=scqx(i,j)*scqy(i,j)
!
scuxi(i,j)=1.0/max(scux(i,j),epsil)
scvyi(i,j)=1.0/max(scvy(i,j),epsil)
scp2i(i,j)=1.0/max(scp2(i,j),epsil)
scq2i(i,j)=1.0/max(scq2(i,j),epsil)
!
! --- largest grid spacing (within limits) used in all diffusion
! --- coefficients: min(max(sc?x,sc?y),sc?x*aspmax,sc?y*aspmax)
aspux(i,j)=min(max(scux(i,j),scuy(i,j)), &
min(scux(i,j),scuy(i,j))*aspmax) &
/max(scux(i,j),epsil)
aspuy(i,j)=min(max(scux(i,j),scuy(i,j)), &
min(scux(i,j),scuy(i,j))*aspmax) &
/max(scuy(i,j),epsil)
aspvx(i,j)=min(max(scvx(i,j),scvy(i,j)), &
min(scvx(i,j),scvy(i,j))*aspmax) &
/max(scvx(i,j),epsil)
aspvy(i,j)=min(max(scvx(i,j),scvy(i,j)), &
min(scvx(i,j),scvy(i,j))*aspmax) &
/max(scvy(i,j),epsil)
!
util1(i,j)=depths(i,j)*scp2(i,j)
enddo
enddo
!
! --- read ice shelf depth array
!
if (ishelf.eq.0) then
ishlf(:,:) = ip(:,:) !no ice shelf
else
if (mnproc.eq.1) then ! .b file from 1st tile only
write (lp,'(3a)') ' reading ice shelf file from ', &
trim(flnmshlf),'.[ab]'
open (unit=uoff+9,file=trim(flnmshlf)//'.b', &
status='old')
read ( uoff+9,'(a79)') preambl
endif
call xcsync(flush_lp)
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
i = index(cline,'=')
read (cline(i+1:),*) hminb,hmaxb
if (mnproc.eq.1) then ! .b file from 1st tile only
close(unit=uoff+9)
write (lp,'(/(1x,a))') preambl,cline
endif
!
call zaiopf(trim(flnmshlf)//'.a','old', 9)
call zaiord(util3,ip,.false., hmina,hmaxa, 9)
call zaiocl(9)
!
if (abs(hmina-hminb).gt.abs(hminb)*1.e-4 .or. &
abs(hmaxa-hmaxb).gt.abs(hmaxb)*1.e-4 ) then
if (mnproc.eq.1) then
write(lp,'(/ a / a,1p3e14.6 / a,1p3e14.6 /)') &
'error - .a and .b files not consistent:', &
'.a,.b min = ',hmina,hminb,hmina-hminb, &
'.a,.b max = ',hmaxa,hmaxb,hmaxa-hmaxb
endif
call xcstop('(geopar)')
stop '(geopar)'
endif
!
!$OMP PARALLEL DO PRIVATE(j,i) &
!$OMP SCHEDULE(STATIC,jblk)
do j= 1,jj
do i= 1,ii
if (ip(i,j).eq.0) then
util3(i,j) = 0.0 !land
elseif (util3(i,j).gt.0.5*hugel) then
util3(i,j) = 0.0 !ice shelf over ocean
elseif (util3(i,j).le.0.0) then
util3(i,j) = 0.0 !ice shelf over ocean
else
util3(i,j) = 1.0 !open ocean
endif
enddo
enddo
call xctilr(util3,1,1, nbdy,nbdy, halo_ps)
ishlf(:,:) = 0 !for jj:jdm and ii:idm
!$OMP PARALLEL DO PRIVATE(j,i) &
!$OMP SCHEDULE(STATIC,jblk)
do j= 1-nbdy,jj+nbdy
do i= 1-nbdy,ii+nbdy
ishlf(i,j) = util3(i,j)
util2(i,j) = ip(i,j)
enddo
enddo
!
call xcsum(sum_is, util3,ip)
call xcsum(sum_ip, util2,ip)
call xcsum(sum_isa, scp2, ishlf)
call xcsum(area, scp2, ip)
nishlf = nint(sum_ip) - nint(sum_is)
if (mnproc.eq.1) then
write (lp,'(/a,i9,f10.2)') &
' number of ice shelf points and area (10^6 km^2):', &
nishlf,(area-sum_isa)*1.d-12
endif
call xcsync(flush_lp)
endif !ishelf
!
! --- In arctic (tripole) domain, top row of mass points is redundent,
! --- so always use ipa, based on ishlf, for mass sums
#if defined(ARCTIC)
ipa(:,:) = ishlf(:,:)
if (jj+j0.eq.jtdm) then
! --- mask top row of mass points
ipa(:,jj:jj+nbdy) = 0
endif
#else
! --- Not a tripole domain, so ipa=ishlf
ipa(:,:) = ishlf(:,:)
#endif
!
call xcsum(avgbot, util1,ipa)
call xcsum(area, scp2, ipa)
avgbot=avgbot/area
if (mnproc.eq.1) then
write (lp,'(/a,f9.1,f10.2)') &
' mean basin depth (m) and area (10^6 km^2):', &
avgbot,area*1.e-12
endif
call xcsync(flush_lp)
!
! --- calculate dx0k?
if (dx00.gt.0.0) then
! --- logarithmic k-dependence of dx0
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'dx00 =',dx00
write(lp,*) 'dx00x =',dx00x
write(lp,*) 'dx00f =',dx00f
endif
call xcsync(flush_lp)
!
dx00 =onem*dx00
dx00x=onem*dx00x
dx0k(1)=dx00
!
dx0kf=1.0
do k=2,kk
dx0kf=dx0kf*dx00f
if (mnproc.eq.1) then
write(lp,*) 'dx0kf =',dx0kf,k
endif
call xcsync(flush_lp)
dx0k(k)=min(dx00*dx0kf,dx00x)
enddo !k
else
do k=1,kk
dx0k(k)=onem*dx0k(k)
enddo !k
endif
if (mnproc.eq.1) then
write(lp,*)
write(lp,125) 1,dx0k( 1)*qonem
write(lp,125) kk,dx0k(kk)*qonem
endif
125 format('dx0k(',i2,') =',f7.2,' m')
!
! --- calculate dp0k and ds0k?
if (dp00.lt.0.0) then
! --- dp0k and ds0k already input
dp00 =onem*dp0k(1)
dp00x=onem*dp0k(kk-1)
dp00i=onem*dp00i
dpms = 0.0
do k=1,kk
dpm = dp0k(k)
dpms = dpms + dpm
dp0k(k) = dp0k(k)*onem
if (mnproc.eq.1) then
write(lp,135) k,dp0k(k)*qonem,dpm,dpms
endif
if (mnproc.eq.-99) then ! bugfix that prevents optimization
write(6,*) 'geopar: dp0k = ',dp0k(k),k,mnproc
endif
call xcsync(flush_lp)
enddo !k
dsms = 0.0
do k=1,nsigma
dsm = ds0k(k)
dsms = dsms + dsm
ds0k(k) = ds0k(k)*onem
if (mnproc.eq.1) then
write(lp,130) k,ds0k(k)*qonem,dsm,dsms
endif
if (mnproc.eq.-99) then ! bugfix that prevents optimization
write(6,*) 'geopar: ds0k = ',ds0k(k),k,mnproc
endif
call xcsync(flush_lp)
enddo !k
if (mnproc.eq.1) then
write(lp,*)
endif
else
! --- calculate dp0k and ds0k
!
! --- logarithmic k-dependence of dp0 (deep z's)
dp00 =onem*dp00
dp00x=onem*dp00x
dp00i=onem*dp00i
if (isopyc) then
dp0k(1)=thkmin*onem
else
dp0k(1)=dp00
endif
dpm = dp0k(1)*qonem
dpms = dpm
if (mnproc.eq.1) then
write(lp,*)
write(lp,135) 1,dp0k(1)*qonem,dpm,dpms
endif
135 format('dp0k(',i2,') =',f7.2,' m', &
' thkns =',f7.2,' m', &
' depth =',f8.2,' m')
call xcsync(flush_lp)
!
dp0kf=1.0
do k=2,kk
dp0kf=dp0kf*dp00f
if (k.le.nhybrd) then
if (dp00f.ge.1.0) then
dp0k(k)=min(dp00*dp0kf,dp00x)
else
dp0k(k)=max(dp00*dp0kf,dp00x)
endif
else
dp0k(k)=0.0
endif
dpm = dp0k(k)*qonem
dpms = dpms + dpm
if (mnproc.eq.1) then
write(lp,135) k,dp0k(k)*qonem,dpm,dpms
endif
if (mnproc.eq.-99) then ! bugfix that prevents optimization
write(6,*) 'geopar: dp0kf = ',dp0kf, mnproc
write(6,*) 'geopar: dp0k = ',dp0k(k),k,mnproc
endif
call xcsync(flush_lp)
enddo !k
!
! --- logarithmic k-dependence of ds0 (shallow z-s)
ds00 =onem*ds00
ds00x=onem*ds00x
if (isopyc) then
ds0k(1)=thkmin*onem
else
ds0k(1)=ds00
endif
dsm = ds0k(1)*qonem
dsms = dsm
if (mnproc.eq.1) then
write(lp,*)
write(lp,130) 1,ds0k(1)*qonem,dsm,dsms
endif
130 format('ds0k(',i2,') =',f7.2,' m', &
' thkns =',f7.2,' m', &
' depth =',f8.2,' m')
call xcsync(flush_lp)
!
ds0kf=1.0
do k=2,nsigma
ds0kf=ds0kf*ds00f
if (ds00f.ge.1.0) then
ds0k(k)=min(ds00*ds0kf,ds00x)
else
ds0k(k)=max(ds00*ds0kf,ds00x)
endif
dsm = ds0k(k)*qonem
dsms = dsms + dsm
if (mnproc.eq.1) then
write(lp,130) k,ds0k(k)*qonem,dsm,dsms
endif
if (mnproc.eq.-99) then ! bugfix that prevents optimization
write(6,*) 'geopar: ds0kf = ',ds0kf, mnproc
write(6,*) 'geopar: ds0k = ',ds0k(k),k,mnproc
endif
call xcsync(flush_lp)
enddo !k
if (mnproc.eq.1) then
write(lp,*)
endif
endif !input:calculate dp0k,ds0k
!
! --- start and stop depths for terrain following coordinate
if (nsigma.eq.0) then
dpns = dp0k(1)
dsns = 0.0
ds0k(1) = dp0k(1)
do k= 2,kk
ds0k(k)=0.0
enddo !k
else
dpns = 0.0
dsns = 0.0
do k=1,nsigma
dpns = dpns + dp0k(k)
dsns = dsns + ds0k(k)
enddo !k
do k= nsigma+1,kk
ds0k(k)=0.0
enddo !k
endif !nsigma
dpns = dpns*qonem !depths is in m
dsns = dsns*qonem !depths is in m
!
if (mnproc.eq.1) then
write(lp,131) nsigma,dpns,dsns
endif
131 format('nsigma = ',i2, &
' deep =',f8.2,' m', &
' shallow =',f8.2,' m' )
call flush(lp)
!
! --- initialize thermobaric reference state arrays.
!
if (kapref.eq.-1) then
if (mnproc.eq.1) then ! .b file from 1st tile only
write (lp,'(3a)') ' reading thermobaric reference file from ', &
trim(flnmforw), 'tbaric.[ab]'
open (unit=uoff+9,file=trim(flnmforw)//'tbaric.b', &
status='old')
read ( uoff+9,'(a79)') preambl
endif
call xcsync(flush_lp)
call zagetc(cline,ios, uoff+9)
if (ios.ne.0) then
if (mnproc.eq.1) then
write(lp,'(/ a,i4,i9 /)') &
'geopar: I/O error from zagetc, iunit,ios = ',uoff+9,ios
endif !1st tile
call xcstop('(geopar)')
stop '(geopar)'
endif
i = index(cline,'=')
read (cline(i+1:),*) hminb,hmaxb
if (mnproc.eq.1) then ! .b file from 1st tile only
close(unit=uoff+9)
write (lp,'(/(1x,a))') preambl,cline
endif
!
! --- input field is between 1.0 and 3.0 and indicates the
! --- relative strength of the two nearest reference states,
! --- e.g. 1.7 is 70% ref2 and 30% ref1
! --- and 2.3 is 70% ref2 and 30% ref3.
!
call zaiopf(trim(flnmforw)//'tbaric.a','old', 9)
call zaiord(util1,ip,.false., hmina,hmaxa, 9)
call zaiocl(9)
!
if (abs(hmina-hminb).gt.abs(hminb)*1.e-4 .or. &
abs(hmaxa-hmaxb).gt.abs(hmaxb)*1.e-4 ) then
if (mnproc.eq.1) then
write(lp,'(/ a / a,1p3e14.6 / a,1p3e14.6 /)') &
'error - .a and .b files not consistent:', &
'.a,.b min = ',hmina,hminb,hmina-hminb, &
'.a,.b max = ',hmaxa,hmaxb,hmaxa-hmaxb
endif
call xcstop('(geopar)')
stop '(geopar)'
endif
!
do j= 1,jj
do i= 1,ii
if (ip(i,j).eq.0) then
util1(i,j) = 1.0 !land
endif
enddo
enddo
!
vland = 1.0
call xctilr(util1, 1,1, nbdy,nbdy, halo_ps)
vland = 0.0
!
! kapi is the 2nd reference state (1st is always 2)
! skap is the scale factor (0.0-1.0) for the 1st reference state
!
! assumes that reference states 1 and 3 are never next to each other.
!
util2(:,:) = -99.0 !otherwise ii(jj)+1 to i(j)dm are NaN
do j= 1,jj
do i= 1,ii
if (max(util1(i, j), &
util1(i-1,j), &
util1(i+1,j), &
util1(i, j-1), &
util1(i, j+1) ).gt.2.0) then
util2(i,j) = 3.0 !kapi
skap(i,j) = 3.0 - util1(i,j)
else
util2(i,j) = 1.0 !kapi
skap(i,j) = util1(i,j) - 1.0
endif
enddo
enddo
vland = 1.0
call xctilr(util2, 1,1, nbdy,nbdy, halo_ps)
call xctilr(skap, 1,1, nbdy,nbdy, halo_ps)
vland = 0.0
!
kapi(:,:) = util2(:,:)
else
skap(:,:) = 1.0 !for diagnostics only
kapi(:,:) = kapref !for diagnostics only
endif !kapref.eq.-1:else
!
! --- initialize some arrays
! --- set depthu,dpu,utotn,pgfx,depthv,dpv,vtotn,pgfy to zero everywhere,
! --- so that they can be used at "lateral neighbors" of u and v points.
! --- similarly for pbot,dp at neighbors of q points.
!
disp_count=0
!
!$OMP PARALLEL DO PRIVATE(j,i,k,ktr) &
!$OMP SCHEDULE(STATIC,jblk)
do j=1-nbdy,jj+nbdy
do i=1-nbdy,ii+nbdy
p( i,j,1)=0.0
pu( i,j,1)=0.0
pv( i,j,1)=0.0
utotn( i,j)=0.0
vtotn( i,j)=0.0
pgfx( i,j)=0.0
pgfy( i,j)=0.0
gradx( i,j)=0.0
grady( i,j)=0.0
depthu(i,j)=0.0
depthv(i,j)=0.0
pbot( i,j)=0.0
!
displd_mn(i,j)=0.0
dispqd_mn(i,j)=0.0
tidepg_mn(i,j)=0.0
!
psikk( i,j,1)=0.0
psikk( i,j,2)=0.0
thkk( i,j,1)=0.0
thkk( i,j,2)=0.0
!
ubavg( i,j,1)=hugel
ubavg( i,j,2)=hugel
ubavg( i,j,3)=hugel
vbavg( i,j,1)=hugel
vbavg( i,j,2)=hugel
vbavg( i,j,3)=hugel
!
umix( i,j)=hugel
vmix( i,j)=hugel
utotm( i,j)=hugel
vtotm( i,j)=hugel
uflux( i,j)=hugel
vflux( i,j)=hugel
uflux2(i,j)=hugel
vflux2(i,j)=hugel
uflux3(i,j)=hugel
vflux3(i,j)=hugel
uja( i,j)=hugel
ujb( i,j)=hugel
via( i,j)=hugel
vib( i,j)=hugel
do k=1,kk
dp( i,j,k,1)=0.0
dp( i,j,k,2)=0.0
dpo(i,j,k,1)=0.0
dpo(i,j,k,2)=0.0
dpu(i,j,k,1)=0.0
dpu(i,j,k,2)=0.0
dpv(i,j,k,1)=0.0
dpv(i,j,k,2)=0.0
!
u( i,j,k,1)=hugel
u( i,j,k,2)=hugel
v( i,j,k,1)=hugel
v( i,j,k,2)=hugel
!
uflx( i,j,k)=hugel
vflx( i,j,k)=hugel
!
dpav( i,j,k)=0.0
uflxav(i,j,k)=0.0
vflxav(i,j,k)=0.0
diaflx(i,j,k)=0.0
!
do ktr= 1,ntracr
tracer(i,j,k,1,ktr)=0.0
tracer(i,j,k,2,ktr)=0.0
enddo
enddo
enddo
enddo
!$OMP END PARALLEL DO
!
!$OMP PARALLEL DO PRIVATE(j,l,i,k) &
!$OMP SCHEDULE(STATIC,jblk)
do j=1,jj
do l=1,isp(j) !ok
do i=max(1,ifp(j,l)),min(ii,ilp(j,l)+1)
ubavg(i,j,1)=0.0
ubavg(i,j,2)=0.0
ubavg(i,j,3)=0.0
umix (i,j)=0.0
utotm (i,j)=0.0
uflux (i,j)=0.0
uflux2(i,j)=0.0
uflux3(i,j)=0.0
uja(i,j)=0.0
ujb(i,j)=0.0
!
do k=1,kk
uflx(i,j,k)=0.0
u(i,j,k,1)=0.0
u(i,j,k,2)=0.0
enddo
enddo
enddo
enddo
!
call xctilr(ubavg, 1, 3, nbdy,nbdy, halo_us) ! note scalar
call xctilr(umix, 1, 1, nbdy,nbdy, halo_us) ! note scalar
call xctilr(utotm, 1, 1, nbdy,nbdy, halo_us) ! note scalar
call xctilr(uflux, 1, 1, nbdy,nbdy, halo_us) ! note scalar
call xctilr(uflux2, 1, 1, nbdy,nbdy, halo_us) ! note scalar
call xctilr(uflux3, 1, 1, nbdy,nbdy, halo_us) ! note scalar
call xctilr(uja, 1, 1, nbdy,nbdy, halo_us)
call xctilr(ujb, 1, 1, nbdy,nbdy, halo_us)
call xctilr(uflx, 1, kk, nbdy,nbdy, halo_us) ! note scalar
call xctilr(u, 1,2*kk, nbdy,nbdy, halo_us) ! note scalar
!
!$OMP PARALLEL DO PRIVATE(i,l,j,k) &
!$OMP SCHEDULE(STATIC)
do i=1,ii
do l=1,jsp(i) !ok
do j=max(1,jfp(i,l)),min(jj,jlp(i,l)+1)
vbavg(i,j,1)=0.0
vbavg(i,j,2)=0.0
vbavg(i,j,3)=0.0
vmix (i,j)=0.0
vtotm (i,j)=0.0
vflux (i,j)=0.0
vflux2(i,j)=0.0
vflux3(i,j)=0.0
via(i,j)=0.0
vib(i,j)=0.0
!
do k=1,kk
vflx(i,j,k)=0.0
v(i,j,k,1)=0.0
v(i,j,k,2)=0.0
enddo
enddo
enddo
enddo
!
call xctilr(vbavg, 1, 3, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vmix, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vtotm, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vflux, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vflux2, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vflux3, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(via, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vib, 1, 1, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(vflx, 1, kk, nbdy,nbdy, halo_vs) ! note scalar
call xctilr(v, 1,2*kk, nbdy,nbdy, halo_vs) ! note scalar
!
return
end
subroutine geopar_halo(sc, g_sc, halo_type)
use mod_xc ! HYCOM communication interface
use mod_cb_arrays ! HYCOM saved arrays
implicit none
!
integer halo_type
real sc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy)
real*4 g_sc(itdm,jtdm)
!
! --- define the halo of sc, from g_sc
! --- needed because xctilr may use vland in the halo
!
integer i,ia,ig,j,ja,jg
!
#if defined(ARCTIC)
! --- periodic in x, closed in South, dipole patch in North
!
do j=1,jj
jg = j0+j
do i=1-nbdy,0
ig = mod(i0+i-1+itdm,itdm)+1
sc(i,j) = g_sc(ig,jg)
enddo !i
do i=ii+1,ii+nbdy
ig = mod(i0+i-1+itdm,itdm)+1
sc(i,j) = g_sc(ig,jg)
enddo !i
enddo !j
!
do i=1-nbdy,ii+nbdy
ig = mod(i0+i-1+itdm,itdm)+1
do j=1-nbdy,0
jg = max(1,j0+j)
sc(i,j) = g_sc(ig,jg)
enddo !j
do j=jj+1,jj+nbdy
jg = j0+j
if (jg.le.jtdm) then
sc(i,j) = g_sc(ig,jg)
else
if (halo_type.eq.halo_ps) then
ja = jtdm-1-(jg-jtdm)
ia = itdm-mod(ig-1,itdm)
sc(i,j) = g_sc(ia,ja)
else !halo_type.eq.halo_qs
ja = jtdm-(jg-jtdm)
ia = mod(itdm-(ig-1),itdm)+1
sc(i,j) = g_sc(ia,ja)
endif
endif
enddo !j
enddo !i
#else
do j=1,jj
jg = j0+j
do i=1-nbdy,0
if (nreg.le.2) then
ig = max( 1,i0+i) !closed
else
ig = mod(i0+i-1+itdm,itdm)+1 !periodic
endif
sc(i,j) = g_sc(ig,jg)
enddo !i
do i=ii+1,ii+nbdy
if (nreg.le.2) then
ig = min(itdm,i0+i) !closed
else
ig = mod(i0+i-1+itdm,itdm)+1 !periodic
endif
sc(i,j) = g_sc(ig,jg)
enddo !i
enddo !j
!
do i=1-nbdy,ii+nbdy
if (nreg.le.2) then
ig = max(1, min( itdm, i0+i ) ) !closed
else
ig = mod(i0+i-1+itdm,itdm)+1 !periodic
endif
do j=1-nbdy,0
if (nreg.eq.0 .or. nreg.eq.4) then
jg = max( 1,j0+j) !closed
else
jg = mod(j0+j-1+jtdm,jtdm)+1 !periodic
endif
sc(i,j) = g_sc(ig,jg)
enddo !j
do j=jj+1,jj+nbdy
if (nreg.eq.0 .or. nreg.eq.4) then
jg = min(jtdm,j0+j) !closed
else
jg = mod(j0+j-1+jtdm,jtdm)+1 !periodic
endif
sc(i,j) = g_sc(ig,jg)
enddo !j
enddo !i
#endif
!diag
!diag if (itest.gt.0 .and. jtest.gt.0) then
!diag do j= 1,jtdm
!diag write(lp,'(a,2i5,e16.8)') &