-
Notifications
You must be signed in to change notification settings - Fork 0
/
Copy pathmod_incupd.F
executable file
·796 lines (789 loc) · 24.2 KB
/
mod_incupd.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
#if defined(ROW_LAND)
#define SEA_P .true.
#define SEA_U .true.
#define SEA_V .true.
#elif defined(ROW_ALLSEA)
#define SEA_P allip(j).or.ip(i,j).ne.0
#define SEA_U alliu(j).or.iu(i,j).ne.0
#define SEA_V alliv(j).or.iv(i,j).ne.0
#else
#define SEA_P ip(i,j).ne.0
#define SEA_U iu(i,j).ne.0
#define SEA_V iv(i,j).ne.0
#endif
module mod_incupd
use mod_xc ! HYCOM communication interface
c
implicit none
c
c --- HYCOM incremental updating (for data assimilation)
c
integer, save, public ::
& incflg, ! incremental update flag (0=no,1=yes,2=full-velocity)
& incstp, ! no. timesteps for full update (1=full insertion)
& incupf, ! number of days of incremental updating input
& incice ! direct insertion of sea ice concentration flag
c
integer, save, private ::
& ncount, ! increment time step counter
& ncountd ! increment day counter
c
real*8, save, private ::
& dtimeu ! next days increment field
c
real, allocatable, dimension(:,:),
& save, private ::
& ubinc, ! ubaro increment
& vbinc ! vbaro increment
c
real, allocatable, dimension(:,:,:),
& save, private ::
& tinc, ! t increment
& sinc, ! s increment
& dpinc, ! dp increment
& uinc, ! u increment
& vinc ! v increment
contains
subroutine incupd_init(dtime0)
c
real*8 dtime0
c
c --- subroutine used to calculate increment field for the incremental updating
c --- version: dec 2005
c
integer i,j,k
logical lopen
c
c --- allocate arrays
c
allocate( tinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy,kdm),
& sinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy,kdm),
& dpinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy,kdm),
& uinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy,kdm),
& vinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy,kdm),
& ubinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy),
& vbinc(1-nbdy:idm+nbdy,1-nbdy:jdm+nbdy) )
call mem_stat_add( 5*(idm+2*nbdy)*(jdm+2*nbdy)*kdm )
call mem_stat_add( 2*(idm+2*nbdy)*(jdm+2*nbdy) )
c
c --- set counter to zero
c
ncount=0
ncountd=0
dtimeu=1.d-6
c
c --- read the target fields, and initialize the "inc" arrays.
c
call incupd_read(dtime0)
c
return
end subroutine incupd_init
subroutine incupd_rd(dtime0)
c
real*8 dtime0
c
c --- subroutine used to calculate increment field for the incremental updating
c --- version: dec 2005
c
integer i,j,k
logical lopen
c
if (ncountd.gt.incupf) then
if (ncountd.eq.incupf+1) then
c should never get here
if (mnproc.eq.1) then
write(lp,*) '... ended updating fields with increments ....'
write(lp,*) 'ncountd= ',ncountd
write(lp,*)
endif !1st tile
call xcsync(flush_lp)
endif
return
endif
c
c --- read the target fields, and initialize the "inc" arrays.
c
call incupd_read(dtime0)
c
return
end subroutine incupd_rd
subroutine incupd(n)
use mod_cb_arrays ! HYCOM saved arrays
use mod_pipe ! HYCOM debugging interface
C
integer n
c
c**********
c*
c 1) update hycom variables with increments.
c
c 2) parameters:
c
c output:
c incremental updated model variables
c
c 4) Ole Martin Smedstad (PSI), December 2005
c
c**********
c
logical, parameter :: lpipe_incupd=.false. !extra checking
c
character utxt*12,vtxt*12
integer i,j,k
real utotij,vtotij
c
include 'stmt_fns.h'
c
if (ncountd.gt.incupf) then
c-------if (mnproc.eq.1) then
c-------write(lp,*) '... ended updating fields with increments .....'
c-------write(lp,*) 'ncountd= ',ncountd
c-------write(lp,*)
c-------endif !1st tile
c-------call xcsync(flush_lp)
return
endif
c
c --- update counter
c
if (incstp.ne.1) then
ncount=ncount+1
endif
c
if (ncount.gt.incstp) then
if (ncount.eq.incstp+1) then
if (mnproc.eq.1) then
write(lp,*) '... ended updating fields with increments ...'
write(lp,*) 'ncount= ',ncount
write(lp,*)
endif !1st tile
call xcsync(flush_lp)
endif !ncount==incstp+1
return
endif !ncount>incstp
c
if (mnproc.eq.1) then
write(lp,*)
if (incflg.eq.1) then
write(lp,'(2a)') 'update fields with increments, ',
& 'but not ubavg and vbavg'
else !incflg.eq.2
write(lp,'(2a)') 'update fields with increments, ',
& 'including ubavg and vbavg'
endif !incflg
write(lp,*) '..........ncount= ',ncount
endif !1st tile
call xcsync(flush_lp)
c
c --- incremental update of dp (dpu, dpv).
c
!$OMP PARALLEL DO PRIVATE(j,k,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
if (SEA_P) then
do k=1,kk-1
c --- dp must be non-negative
dp(i,j,k,n) = max( dp(i,j,k,n) + dpinc(i,j,k), 0.0 )
c --- p must be at or above the bottom
p(i,j,k+1) = min( p(i,j,k) + dp(i,j,k,n),
& p(i,j,kk+1) )
dp(i,j,k,n) = p(i,j,k+1) - p(i,j,k)
enddo !k
c --- layer kk always touches the bottom
dp(i,j,kk,n) = p(i,j,kk+1) - p(i,j,kk)
endif !ip
enddo !i
enddo !j
!$OMP END PARALLEL DO
c
call dpudpv(dpu(1-nbdy,1-nbdy,1,n),
& dpv(1-nbdy,1-nbdy,1,n),
& p,depthu,depthv, 0,0)
c
if (lpipe_incupd) then
do k= 1,kk
write (utxt,'(a9,i3)') 'up dpu k=',k
write (vtxt,'(a9,i3)') 'up dpv k=',k
call pipe_compare_sym2(dpu(1-nbdy,1-nbdy,1,n),iu,utxt,
& dpv(1-nbdy,1-nbdy,1,n),iv,vtxt)
enddo !k
endif !lpipe_incupd
c
c --- incremental update of the other fields.
c --- salinity from updated th&S.
c --- rebalance u and v via utotij and vtotij.
c
!$OMP PARALLEL DO PRIVATE(j,i,k,utotij,vtotij)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
if (SEA_P) then
do k=1,kk
if (tinc(i,j,k).ne.0.0 .or.
& sinc(i,j,k).ne.0.0 ) then
temp(i,j,k,n) = temp(i,j,k,n) + tinc(i,j,k)
saln(i,j,k,n) = saln(i,j,k,n) + sinc(i,j,k)
th3d(i,j,k,n) = sig(temp(i,j,k,n),saln(i,j,k,n))-thbase
endif !non-zero increment
enddo ! k
endif !ip
if (SEA_U) then
utotij = 0.0
do k=1,kk
u(i,j,k,n) = u(i,j,k,n) + uinc(i,j,k)
utotij = utotij + u(i,j,k,n)*dpu(i,j,k,n)
enddo ! k
utotij=utotij/depthu(i,j)
do k=1,kk
u(i,j,k,n) = u(i,j,k,n) - utotij
enddo ! k
if (incflg.eq.2) then !update ubavg
ubavg(i,j,n) = ubavg(i,j,n) + ubinc(i,j)
* ubavg(i,j,n) = ubavg(i,j,n) + ubinc(i,j) + utotij
endif !incflg==2
endif !iu
if (SEA_V) then
vtotij = 0.0
do k=1,kk
v(i,j,k,n) = v(i,j,k,n) + vinc(i,j,k)
vtotij = vtotij + v(i,j,k,n)*dpv(i,j,k,n)
enddo ! k
vtotij=vtotij/depthv(i,j)
do k=1,kk
v(i,j,k,n) = v(i,j,k,n) - vtotij
enddo ! k
if (incflg.eq.2) then !update vbavg
vbavg(i,j,n) = vbavg(i,j,n) + vbinc(i,j)
* vbavg(i,j,n) = vbavg(i,j,n) + vbinc(i,j) + vtotij
endif !incflg==2
endif !iv
enddo !i
enddo ! j
!$OMP END PARALLEL DO
c
if (mnproc.eq.1) then
write(lp,*) 'finished incupdate',ncount
write(lp,*)
endif !1st tile
call xcsync(flush_lp)
c
return
end subroutine incupd
subroutine incupd_read(dtime)
use mod_cb_arrays ! HYCOM saved arrays
use mod_za ! HYCOM I/O interface
use mod_pipe ! HYCOM debugging interface
c
real*8 dtime
c
c --- input 3-d HYCOM fields (from an archive file) on model day dtime.
c --- directly insert the input covice and thkice (if they exist).
c --- calculate the increment between the input and the initial state.
c
c --- filenames incup/incupd.iyear_iday_ihour.[ab].
c --- I/O and array I/O unit 925 used here, but not reserved.
c
logical, parameter :: ldebug_incupd_read=.false. !usually .false.
c
character flnm*24, cline*80, cvarin*6, cfield*8
character ptxt*12,utxt*12,vtxt*12
integer i,idmtst,ios,j,jdmtst,k,l,layer,nskip
integer iyear,iday,ihour
logical nodens
real tincstp
* real sumdpi
c
integer nstep0
real*8 dtime0
c
include 'stmt_fns.h'
c
call forday(dtime, yrflag, iyear,iday,ihour)
c
write(flnm,'("incup/incupd.",i4.4,"_",i3.3,"_",i2.2)')
& iyear,iday,ihour
c
if(dtime.ge.dtimeu) then
c
ncountd=ncountd+1
ncount=0
c
if (ncountd.gt.incupf) then
if (ncountd.eq.incupf+1) then
if (mnproc.eq.1) then
write(lp,*) '... ended updating fields with increments ...'
write(lp,*) 'ncountd= ',ncountd
write(lp,*)
endif !1st tile
call xcsync(flush_lp)
ncountd=incupf+99 !turn off "ended" printout
endif
return
endif !ncountd>incupf
c
if (mnproc.eq.1) then
write(lp,*) 'read incremental updating ...'
write(lp,*) 'ncountd ...',ncountd
write (lp,*) 'incupd_read: ',flnm
write (lp,*) ' time: ',dtime
write (lp,*) 'iyear,iday,ihour: ',iyear,iday,ihour
endif !1st tile
call xcsync(flush_lp)
c
call zaiopf(flnm//'.a','old', 925)
if (mnproc.eq.1) then ! .b file from 1st tile only
open (unit=uoff+925,file=flnm//'.b',form='formatted',
& status='old',action='read')
c
read(uoff+925,'(a)') cline
read(uoff+925,'(a)') cline
read(uoff+925,'(a)') cline
read(uoff+925,'(a)') cline
c
read(uoff+925,'(a)') cline
read(uoff+925,'(a)') cline
read(uoff+925,'(a)') cline
endif !1st tile
c
call zagetc(cline,ios, uoff+925)
read(cline,*) idmtst,cvarin
* if (mnproc.eq.1) then
* write(lp,*) cvarin,' = ',idmtst
* endif !1st tile
if (cvarin.ne.'idm ') then
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'error in incupd_read - input ',cvarin,
& ' but should be idm '
write(lp,*)
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
call zagetc(cline,ios, uoff+925)
read(cline,*) jdmtst,cvarin
* if (mnproc.eq.1) then
* write(lp,*) cvarin,' = ',jdmtst
* endif !1st tile
if (cvarin.ne.'jdm ') then
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'error in incupd_read - input ',cvarin,
& ' but should be jdm '
write(lp,*)
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
c
if (idmtst.ne.itdm .or. jdmtst.ne.jtdm) then
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'error in incupd_read - input idm,jdm',
& ' not consistent with parameters'
write(lp,*) 'idm,jdm = ',itdm, jtdm, ' (dimensions.h)'
write(lp,*) 'idm,jdm = ',idmtst,jdmtst,' (input)'
write(lp,*)
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
c
if (mnproc.eq.1) then ! .b file from 1st tile only
read (uoff+925,*)
endif
c
c --- skip (most) surface fields.
c
call zaiosk(925)
call zagetc(cline,ios, uoff+925)
i = index(cline,'=')
read(cline(i+1:),*) nstep0,dtime0,layer
if (mnproc.eq.1) then
write(lp,*) 'dtime0= ',dtime0
endif
if (dtime0.ne.dtime) then
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'error in incupd_read - input ',dtime0,
& ' but dtime should be ',dtime
write(lp,*)
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
nodens = layer.ne.0 !new or original archive type
if (nodens .and. layer.ne.sigver) then
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'error in incupd_read - input ',layer,
& ' sigver but should be ',sigver
write(lp,*)
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
c
c assumes that there is a new incremental updating file once a day
c for "incupf" days, see blkdat.input
c
dtimeu=dtime0+1.d0
c
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'dtime, dtime0, dtimeu = ',dtime,
& dtime0, dtimeu
write(lp,*)
endif !1st tile
call xcsync(flush_lp)
c
if (nodens) then
do i= 2,6
if (mnproc.eq.1) then ! .b file from 1st tile only
read (uoff+925,*)
endif
call zaiosk(925)
enddo
else
do i= 2,11
if (mnproc.eq.1) then ! .b file from 1st tile only
read (uoff+925,*)
endif
call zaiosk(925)
enddo
endif
c
call rd_archive(ubinc, cfield,layer, 925) !u_btrop or covice or mix_dpth
if (cfield.eq.'mix_dpth') then
c --- archive contains 'steric '
call rd_archive(ubinc, cfield,layer, 925) !u_btrop or covice
endif
if (mnproc.eq.1) then
write(lp,'(2a)') "surface: ",cfield
endif
call xcsync(flush_lp)
if (cfield.eq.'covice ') then
c
c --- directly insert covice and thkice.
c
call rd_archive(util5, cfield,layer, 925) !thkice
if (mnproc.eq.1) then
write(lp,'(2a)') "surface: ",cfield
endif
call xcsync(flush_lp)
!$OMP PARALLEL DO PRIVATE(j,k,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
if (SEA_P) then
covice(i,j)=ubinc(i,j)
thkice(i,j)=util5(i,j)
endif !ip
enddo !li
enddo !j
!$OMP END PARALLEL DO
call zaiosk(925) !temice
if (mnproc.eq.1) then ! .b file from 1st tile only
read (uoff+925,*)
endif
call rd_archive(ubinc, cfield,layer, 925)
if (mnproc.eq.1) then
write(lp,'(2a)') "surface: ",cfield
endif
call xcsync(flush_lp)
incice = 1 !have input covice, don't direct insert si_c
else
incice = -1 !have not input covice, might direct insert si_c
endif
call rd_archive(vbinc, cfield,layer, 925)
if (mnproc.eq.1) then
write(lp,'(2a)') "surface: ",cfield
endif
call xcsync(flush_lp)
c
if (mnproc.eq.1) then
write (lp,*) 'start 3-D archive file read'
endif
call xcsync(flush_lp)
c
c --- 3-d fields.
c
nskip = 0
do k=1,kk
call rd_archive(uinc(1-nbdy,1-nbdy,k), cfield,layer, 925)
if (cfield.ne.'u-vel. ' .and. k.ne.2) then
if (mnproc.eq.1) then
write(lp,'(/ a / a,a /)') cfield,
& 'error in incupd_read - expected ','u-vel. '
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
elseif (cfield.ne.'u-vel. ') then !k==2
c
c --- count "tracer" fields (to be skipped)
c
if (mnproc.eq.1) then
write(lp,'(2a)') "counting tracers: ",cfield
endif
do nskip= 2,99
call rd_archive(uinc(1-nbdy,1-nbdy,k), cfield,layer, 925)
if (mnproc.eq.1) then
write(lp,'(2a)') "counting tracers: ",cfield
endif
if (cfield.eq.'u-vel. ') then
exit
endif
enddo !nskip
nskip = nskip - 1
write(lp,'(a,i3)') "nskip =",nskip
endif
call rd_archive(vinc(1-nbdy,1-nbdy,k), cfield,layer, 925)
if (cfield.ne.'v-vel. ') then
if (mnproc.eq.1) then
write(lp,'(/ a / a,a /)') cfield,
& 'error in incupd_read - expected ','v-vel. '
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
c if (mnproc.eq.1) then
c write (lp,*) 'read v-vel archive file'
c endif
call rd_archive(dpinc(1-nbdy,1-nbdy,k), cfield,layer, 925)
if (cfield.ne.'thknss ') then
if (mnproc.eq.1) then
write(lp,'(/ a / a,a /)') cfield,
& 'error in incupd_read - expected ','thknss '
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
c if (mnproc.eq.1) then
c write (lp,*) 'read dpinc archive file'
c endif
call rd_archive(tinc(1-nbdy,1-nbdy,k), cfield,layer, 925)
if (cfield.ne.'temp ') then
if (mnproc.eq.1) then
write(lp,'(/ a / a,a /)') cfield,
& 'error in incupd_read - expected ','temp '
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
call rd_archive(sinc(1-nbdy,1-nbdy,k), cfield,layer, 925)
if (cfield.ne.'salin ') then
if (mnproc.eq.1) then
write(lp,'(/ a / a,a /)') cfield,
& 'error in incupd_read - expected ','salin '
endif !1st tile
call xcstop('(incupd_read)')
stop '(incupd_read)'
endif
if (.not. nodens) then
c --- skip density
if (mnproc.eq.1) then ! .b file from 1st tile only
read (uoff+925,*)
endif
call zaiosk(925)
endif !nodens:else
c
c --- skip (nskip) tracers
c
do l= 1,nskip
if (mnproc.eq.1) then ! .b file from 1st tile only
read (uoff+925,*)
endif
call zaiosk(925)
enddo !l
enddo !k
c
call xctilr(dpinc,1,kk, 1,1, halo_ps)
c
if (mnproc.eq.1) then ! .b file from 1st tile only
close( unit=uoff+925)
endif
call zaiocl(925)
c
c --- calculate increments
c --- the "inc" reads, above, are full HYCOM fields (not increments yet).
c
if(incstp.eq.1) then
tincstp=1.0
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'tincstp=1.0 ',tincstp,incstp
endif
else
tincstp=2.0/real(incstp)
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'tincstp=2.0/incstp ',tincstp,incstp
endif
endif !incstp
c
c
if (mnproc.eq.1) then
write(lp,*)
write(lp,*) 'calculate t,s,u,v and dp increments'
endif !1st tile
call xcsync(flush_lp)
c
* if (iutest.gt.0 .and. jutest.gt.0) then
* write(lp,*) '*',' iutest= ',iutest+i0,' jutest= ',jutest+j0,' *'
* write(lp,*) '*********** dpinc input ************'
* sumdpi=0.0
* write(lp,'(a)')
* & 'k,dp1,dp2,dpinc='
* do k= 1,kk
* sumdpi=sumdpi+dpinc(iutest,jutest,k)
* write(lp,'(a,i3,3f20.5)')
* & 'k= ',
* & k,dp(iutest,jutest,k,1)*qonem,
* & dp(iutest,jutest,k,2)*qonem,
* & dpinc(iutest,jutest,k)*qonem
* call flush(lp)
* enddo !k
* write(lp,*) 'sumdpi= ', sumdpi*qonem
* call flush(lp)
* endif
c
!$OMP PARALLEL DO PRIVATE(j,k,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1,jj
do i=1,ii
if (SEA_U) then
ubinc(i,j)=(ubinc(i,j) - ubavg(i,j,1))*tincstp
do k=1,kk
c use an approximate 2*dpu
if (dpinc(i,j,k)+dpinc(i-1,j,k).gt.2.0*onemm) then
uinc(i,j,k)=(uinc(i,j,k) - u(i,j,k,1))*tincstp
else
uinc(i,j,k)=0.0 !thin target layer
endif
enddo !k
endif !iu
if (SEA_V) then
vbinc(i,j)=(vbinc(i,j) - vbavg(i,j,1))*tincstp
do k=1,kk
c use an approximate 2*dpv
if (dpinc(i,j,k)+dpinc(i,j-1,k).gt.2.0*onemm) then
vinc(i,j,k)=(vinc(i,j,k) - v(i,j,k,1))*tincstp
else
vinc(i,j,k)=0.0 !thin target layer
endif
enddo !k
endif !iv
if (SEA_P) then
do k=1,kk
if (dpinc(i,j,k).gt.onemm) then
sinc(i,j,k)=(sinc(i,j,k) - saln(i,j,k,1))*tincstp
tinc(i,j,k)=(tinc(i,j,k) - temp(i,j,k,1))*tincstp
else
tinc(i,j,k)=0.0 !thin target layer
sinc(i,j,k)=0.0 !thin target layer
endif
dpinc(i,j,k)=(dpinc(i,j,k) - dp(i,j,k,1))*tincstp
enddo !k
endif !ip
enddo !li
enddo !j
!$OMP END PARALLEL DO
c
call xctilr(dpinc,1,kk, 1,1, halo_ps)
c
if (ldebug_incupd_read) then
call pipe_compare_sym2(ubinc,iu,'incupd:ubinc',
& vbinc,iv,'incupd:vbinc')
do k= 1,kk
write (utxt,'(a9,i3)') ' uinc k=',k
write (vtxt,'(a9,i3)') ' vinc k=',k
call pipe_compare_sym2(uinc(1-nbdy,1-nbdy,k),iu,utxt,
& vinc(1-nbdy,1-nbdy,k),iv,vtxt)
write (ptxt,'(a9,i3)') ' dpinc k=',k
call pipe_compare_sym1(dpinc(1-nbdy,1-nbdy,k),ip,ptxt)
write (ptxt,'(a9,i3)') ' tinc k=',k
call pipe_compare_sym1( tinc(1-nbdy,1-nbdy,k),ip,ptxt)
write (ptxt,'(a9,i3)') ' sinc k=',k
call pipe_compare_sym1( sinc(1-nbdy,1-nbdy,k),ip,ptxt)
enddo !k
endif !ldebug_incupd_read
*
* if (iutest.gt.0 .and. jutest.gt.0) then
* write(lp,*) '*',' iutest= ',iutest+i0,' jutest= ',jutest+j0,' *'
* write(lp,*) '*********** dpinc out ************'
* write(lp,'(a)')
* & 'k,dp1,dp2,dpinc='
* sumdpi=0.0
* do k= 1,kk
* sumdpi=sumdpi+dpinc(iutest,jutest,k)
* write(lp,'(a,i3,3f20.5)')
* & 'k= ',
* & k,dp(iutest,jutest,k,1)*qonem,
* & dp(iutest,jutest,k,2)*qonem,
* & dpinc(iutest,jutest,k)*qonem
* call flush(lp)
* enddo !k
* write(lp,*) 'inc sumdpi= ', sumdpi*qonem
* call flush(lp)
* endif
c
if (mnproc.eq.1) then
write(lp,*) '... finished reading incupd',dtime,dtime0
endif !1st tile
call xcsync(flush_lp)
c
endif ! dtime
c
return
end subroutine incupd_read
subroutine incupd_si_c(dtime)
use mod_cb_arrays ! HYCOM saved arrays
use mod_za ! HYCOM I/O interface
use mod_pipe ! HYCOM debugging interface
c
real*8 dtime
c
c --- directly insert si_c into covice and thkice.
c
integer i,j,l
c
if (incice.eq.-1) then
incice = 0 !will have directly inserted si_c
c
c --- directly insert covice and thkice.
c
!$OMP PARALLEL DO PRIVATE(j,i)
!$OMP& SCHEDULE(STATIC,jblk)
do j=1-nbdy,jj+nbdy
do i=1,ii
if (SEA_P) then
covice(i,j)=si_c(i,j)
thkice(i,j)=covice(i,j)*hicemn
endif !ip
enddo !i
enddo !j
!$OMP END PARALLEL DO
c
if (mnproc.eq.1) then
write(lp,*) '... finished inserting si_c',dtime
endif !1st tile
call xcsync(flush_lp)
endif !incice
c
return
end subroutine incupd_si_c
c
end module mod_incupd
c
c
c> Revision history:
c>
c> Feb 2006 - 1st module version
c> May 2006 - changed to read multiple increment files
c> Jul 2011 - thin layer is now 1mm (no longer 1m)
c> Jul 2011 - replace thinc with sinc
c> Nov 2012 - bugfix: added xctilr(dpinc to update halo
c> Dec 2012 - cleaned up printing to .log file
c> Apr 2012 - added incice and incupd_si_c
c> May 2014 - use land/sea masks (e.g. ip) to skip land