My Project
mod_obcs2.f90
Go to the documentation of this file.
1 
2 
3 
4 
5 
6 
7 
8 
9 
10 
11 
12 !/===========================================================================/
13 ! Copyright (c) 2007, The University of Massachusetts Dartmouth
14 ! Produced at the School of Marine Science & Technology
15 ! Marine Ecosystem Dynamics Modeling group
16 ! All rights reserved.
17 !
18 ! FVCOM has been developed by the joint UMASSD-WHOI research team. For
19 ! details of authorship and attribution of credit please see the FVCOM
20 ! technical manual or contact the MEDM group.
21 !
22 !
23 ! This file is part of FVCOM. For details, see http://fvcom.smast.umassd.edu
24 ! The full copyright notice is contained in the file COPYRIGHT located in the
25 ! root directory of the FVCOM code. This original header must be maintained
26 ! in all distributed versions.
27 !
28 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS"
29 ! AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO,
30 ! THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
31 ! PURPOSE ARE DISCLAIMED.
32 !
33 !/---------------------------------------------------------------------------/
34 ! CVS VERSION INFORMATION
35 ! $Id$
36 ! $Name$
37 ! $Revision$
38 !/===========================================================================/
39 
40 MODULE mod_obcs2
41 
42  USE all_vars
43  USE mod_prec
44  USE mod_obcs
45  USE mod_types
46  USE bcs
47 
48  USE mod_meanflow
49 
50  IMPLICIT NONE
51  SAVE
52 
53  INTEGER :: nobe,nobcv
54  INTEGER,ALLOCATABLE :: nobedge_lst(:),cobedge_lst(:)
55  INTEGER,ALLOCATABLE :: i_obc_cell(:),i_obc_node(:),i_obc_cell2(:),i_obc_node2(:)
56  REAL(sp),ALLOCATABLE :: uatts(:,:),vatts(:,:),utts(:,:,:),vtts(:,:,:),eltts(:,:)
57 
58  REAL(sp),ALLOCATABLE :: elt(:), eltf(:), elrkt(:), eltdt(:)
59  REAL(sp),ALLOCATABLE :: elp(:), elpf(:), elrkp(:)
60  REAL(sp),ALLOCATABLE :: uat(:), vat(:), uatf(:), vatf(:)
61  REAL(sp),ALLOCATABLE :: uap(:), vap(:)
62 
63  REAL(sp),ALLOCATABLE :: uant (:), vant (:), uan (:), van (:)
64  REAL(sp),ALLOCATABLE :: uanp (:), vanp (:)
65 
66  REAL(sp),ALLOCATABLE :: ut(:,:), vt(:,:)
67  REAL(sp),ALLOCATABLE :: unt(:,:), vnt(:,:), un(:,:), vn(:,:)
68 
69  REAL(sp),ALLOCATABLE :: uapf (:), vapf (:)
70  REAL(sp),ALLOCATABLE :: uantf (:), vantf (:), uanf (:), vanf (:)
71  REAL(sp),ALLOCATABLE :: uanpf (:), vanpf (:)
72  REAL(sp),ALLOCATABLE :: uarknt(:), varknt(:), uarkn(:), varkn(:)
73  REAL(sp),ALLOCATABLE :: untb(:,:), vntb(:,:), unb(:,:), vnb(:,:)
74 
77  INTEGER, ALLOCATABLE :: i_tidenode_gl(:),i_tidenode_n(:)
78  INTEGER, ALLOCATABLE :: i_tidecell_gl(:),i_tidecell_n(:)
79 
80  REAL(sp),ALLOCATABLE :: dltn(:)
81  INTEGER :: ios
82  CONTAINS
83 
84 !=========================================================================|
85  SUBROUTINE alloc_obc2_data
86 
87  IMPLICIT NONE
88 
89  ALLOCATE(elt(0:ntidenode)); elt = zero
90  ALLOCATE(eltf(0:ntidenode)); eltf = zero
91  ALLOCATE(elrkt(0:ntidenode)); elrkt = zero
92  ALLOCATE(eltdt(0:ntidenode)); eltdt = zero
93  ALLOCATE(elp(0:ntidenode)); elp = zero
94  ALLOCATE(elpf(0:ntidenode)); elpf = zero
95  ALLOCATE(elrkp(0:ntidenode)); elrkp = zero
96  ALLOCATE(uat(0:ntidecell)); uat = zero
97  ALLOCATE(vat(0:ntidecell)); vat = zero
98  ALLOCATE(uatf(0:ntidecell)); uatf = zero
99  ALLOCATE(vatf(0:ntidecell)); vatf = zero
100  ALLOCATE(uap(0:ntidecell)); uap = zero
101  ALLOCATE(vap(0:ntidecell)); vap = zero
102  ALLOCATE(uant(0: nmfcell)); uant = zero
103  ALLOCATE(vant(0: nmfcell)); vant = zero
104  ALLOCATE(uan(0: nmfcell)); uan = zero
105  ALLOCATE(van(0: nmfcell)); van = zero
106  ALLOCATE(uanp(0: nmfcell)); uanp = zero
107  ALLOCATE(vanp(0: nmfcell)); vanp = zero
108  ALLOCATE(ut(0:ntidecell,1:kbm1)); ut = zero
109  ALLOCATE(vt(0:ntidecell,1:kbm1)); vt = zero
110  ALLOCATE(unt(0: nmfcell,1:kbm1)); unt = zero
111  ALLOCATE(vnt(0: nmfcell,1:kbm1)); vnt = zero
112  ALLOCATE(un(0: nmfcell,1:kbm1)); un = zero
113  ALLOCATE(vn(0: nmfcell,1:kbm1)); vn = zero
114 
115  ALLOCATE(uapf(0:ntidecell)); uapf = zero
116  ALLOCATE(vapf(0:ntidecell)); vapf = zero
117  ALLOCATE(uantf(0: nmfcell)); uantf = zero
118  ALLOCATE(vantf(0: nmfcell)); vantf = zero
119  ALLOCATE(uanf(0: nmfcell)); uanf = zero
120  ALLOCATE(vanf(0: nmfcell)); vanf = zero
121  ALLOCATE(uanpf(0: nmfcell)); uanpf = zero
122  ALLOCATE(vanpf(0: nmfcell)); vanpf = zero
123  ALLOCATE(uarknt(0: nmfcell)); uarknt = zero
124  ALLOCATE(varknt(0: nmfcell)); varknt = zero
125  ALLOCATE(uarkn(0: nmfcell)); uarkn = zero
126  ALLOCATE(varkn(0: nmfcell)); varkn = zero
127  ALLOCATE(untb(0: nmfcell,1:kbm1)); untb = zero
128  ALLOCATE(vntb(0: nmfcell,1:kbm1)); vntb = zero
129  ALLOCATE(unb(0: nmfcell,1:kbm1)); unb = zero
130  ALLOCATE(vnb(0: nmfcell,1:kbm1)); vnb = zero
131 
132  RETURN
133  END SUBROUTINE alloc_obc2_data
134 !==========================================================================|
135 
136 !========================================================================
137  SUBROUTINE find_obside
139  IMPLICIT NONE
140  INTEGER :: I,ITMP,J,J1,I1,IERR
141  INTEGER :: IA,IB
142  INTEGER,ALLOCATABLE :: TEMP(:)
143  INTEGER,ALLOCATABLE :: NODE_OB(:),CELL_OB(:)
144  INTEGER :: k,NCNT,itemp
145  INTEGER, ALLOCATABLE :: temp1(:)
146  REAL(SP),ALLOCATABLE :: RTEMP1(:,:),RTEMP2(:,:),RTEMP3(:,:,:),RTEMP4(:,:,:)
147  INTEGER :: i2,i3,ii,JN
148  REAL(SP):: DELTX,DELTY,XTMP1,YTMP1,AA1,BB1,CC1,AA2,BB2,CC2
149  REAL(SP)::TTIME
150  ALLOCATE(node_ob(0:mt)); node_ob = 0
151  ALLOCATE(cell_ob(0:nt)); cell_ob = 0
152 !-----------------------------------Jianzhong-------------------------!
153 ! DO I1=1,IOBCN
154 ! J=I_OBC_N(I1)
155 ! J1=NEXT_OBC(I1)
156 ! NODE_OB(J) = 1
157 ! DO I=1,NTVE(J)
158 ! CELL_OB(NBVE(J,I)) = 1
159 ! END DO
160 ! END DO
161  DO i1=1,ibcn(2)
162  jn = obc_lst(2,i1)
163  j=i_obc_n(jn)
164  node_ob(j) = 1
165  DO i=1,ntve(j)
166  cell_ob(nbve(j,i)) = 1
167  END DO
168  END DO
169 !---------------------------------------------------------------------!
170 
171  ALLOCATE(temp(ne)); temp = zero
172  nobe = 0
173 
174  DO i=1,ne
175  ia = iec(i,1)
176  ib = iec(i,2)
177  IF(cell_ob(ia) == 1 .OR. cell_ob(ib) == 1)THEN
178  nobe = nobe + 1
179  temp(nobe) = i
180  END IF
181  END DO
182 
183  ALLOCATE(cobedge_lst(nobe))
184  cobedge_lst(1:nobe) = temp(1:nobe)
185  DEALLOCATE(temp)
186 
187  ALLOCATE(temp(ncv)); temp = zero
188  nobcv = 0
189 
190  DO i=1,ncv
191  ia = niec(i,1)
192  ib = niec(i,2)
193  IF(node_ob(ia) == 1 .OR. node_ob(ib) == 1)THEN
194  nobcv = nobcv + 1
195  temp(nobcv) = i
196  END IF
197  END DO
198 
199  ALLOCATE(nobedge_lst(nobcv))
200  nobedge_lst(1:nobcv) = temp(1:nobcv)
201  DEALLOCATE(temp)
202 
203  DEALLOCATE(node_ob)
204  DEALLOCATE(cell_ob)
205 
206  inmf =45
207  intcell =46
208  intnode =47
209  intelel =48
210  intuv =49
211  CALL fopen(inmf, trim(input_dir)//trim(casename)//'_meanflow.dat' ,"cfr")
212  CALL fopen(intcell,trim(input_dir)//trim(casename)//'_tide_cell.dat' ,"cfr")
213  CALL fopen(intnode,trim(input_dir)//trim(casename)//'_tide_node.dat' ,"cfr")
214  CALL fopen(intelel,trim(input_dir)//trim(casename)//'_tide_el.dat' ,"cfr")
215  CALL fopen(intuv, trim(input_dir)//trim(casename)//'_tide_uv.dat' ,"cfr")
216 
217 !-------------------------------------------------------------------
218 !
219 !----Read in Tidal Current Time Series Data----------------------------------
220 !
221 
222  rewind(intcell)
223  READ(intcell,*) ntidecell_gl
224 
225  ntidecell = 0
226  ntidecell_i = 0
227  IF (ntidecell_gl > 0) THEN
228 
229  ALLOCATE(i_tidecell_gl(ntidecell_gl))
230  DO i=1,ntidecell_gl
231  READ(intcell,*)i_tidecell_gl(i)
232  ENDDO
233  CLOSE(intcell)
234 
235 ! IF(ntidecell_GL > 300) THEN
236 ! WRITE(*,*)'CHANGE FORMAT STATEMENT BELOW TO ACCOMODATE ntidecell_GL='
237 ! WRITE(*,*)ntidecell_GL,' NUMBER OF TIDAL OPEN BOUNDARY CELLS AND RECOMPILE'
238 ! CALL PSTOP
239 ! END IF
240 
241  rewind(intuv)
242 
243 !------------------determine the julian forcing counting-------------------------
244  IF(msr)THEN
245  CALL fopen(111,trim(input_dir)//trim(casename)//'_elj_obc.dat',"cfr")
246  ncnt = 0
247  READ(111,*)
248  READ(111,*)
249  DO WHILE(.true.)
250  READ(111,*,iostat=ios)
251  IF(ios < 0)EXIT
252  ncnt = ncnt + 1
253  END DO
254  IF(ncnt == 0)CALL fatal_error("JULIAN TIDE SELECTED BUT NO DATA IN FILE")
255  END IF
256  elo_tm%NTIMES = ncnt
257  ALLOCATE(elo_tm%TIMES(ncnt))
258  ttime = 0.0_sp
259  DO i=1,ncnt
260  elo_tm%TIMES(i) = ttime
261  ttime = ttime + 720 ! TTIME+DELTT
262  END DO
263 !--------------------------------------------------------------------------------
264  ALLOCATE(rtemp1(ntidecell_gl,elo_tm%NTIMES)) ! assumig we've known ELO_TM
265  ALLOCATE(rtemp2(ntidecell_gl,elo_tm%NTIMES))
266  ALLOCATE(rtemp3(ntidecell_gl,kbm1,elo_tm%NTIMES))
267  ALLOCATE(rtemp4(ntidecell_gl,kbm1,elo_tm%NTIMES))
268  rtemp1 = 0.0_sp
269  rtemp2 = 0.0_sp
270  rtemp3 = 0.0_sp
271  rtemp4 = 0.0_sp
272 
273  IF(msr)THEN
274  DO i=1,elo_tm%NTIMES
275  READ(intuv,'(I7,300f8.4)') itemp,(rtemp1(j,i),j=1,ntidecell_gl)
276  READ(intuv,'(I7,300f8.4)') itemp,(rtemp2(j,i),j=1,ntidecell_gl)
277  DO k = 1,kbm1
278  READ(intuv,'(I7,300f8.4)') itemp,(rtemp3(j,k,i),j=1,ntidecell_gl)
279  READ(intuv,'(I7,300f8.4)') itemp,(rtemp4(j,k,i),j=1,ntidecell_gl)
280  ENDDO
281  END DO
282  END IF
283 
284  close(intuv)
285 
286 !
287 !---Map to Local Domain----------------------------------------
288 !
289  IF(serial) THEN
292  ALLOCATE(i_tidecell_n(ntidecell))
294  ALLOCATE(uatts(ntidecell,elo_tm%NTIMES))
295  ALLOCATE(vatts(ntidecell,elo_tm%NTIMES))
296  ALLOCATE(utts(ntidecell,kbm1,elo_tm%NTIMES))
297  ALLOCATE(vtts(ntidecell,kbm1,elo_tm%NTIMES))
298  uatts = rtemp1
299  vatts = rtemp2
300  utts = rtemp3
301  vtts = rtemp4
302  ENDIF
303 
304 
305  DEALLOCATE(rtemp1,rtemp2,rtemp3,rtemp4)
306 
307 
308  ELSE ! if statement end for ntidecell_GL > 0
309  close(intcell)
310  END IF
311 !-------------------------------------------------------------------
312 !
313 !----Read in Tidal Elevation Time Series Data----------------------------------
314 !
315 
316  rewind(intnode)
317  READ(intnode,*) ntidenode_gl
318 
319  ntidenode = 0
320  ntidenode_i = 0
321  IF (ntidenode_gl > 0) THEN
322 
323  ALLOCATE(i_tidenode_gl(ntidenode_gl))
324  DO i=1,ntidenode_gl
325  READ(intnode,*)i_tidenode_gl(i)
326  ENDDO
327  CLOSE(intnode)
328 
329 ! IF(ntidenode_GL > 300 ) THEN
330 ! WRITE(*,*)'CHANGE FORMAT STATEMENT BELOW TO ACCOMODATE ntidenode_GL=',ntidenode_GL
331 ! WRITE(*,*)' NUMBER OF TIDAL OPEN BOUNDARY NODES AND RECOMPILE'
332 ! CALL PSTOP
333 ! END IF
334 
335  ALLOCATE(rtemp1(ntidenode_gl, elo_tm%NTIMES))
336  rtemp1 = 0.0_sp
337 
338  rewind(intelel)
339  IF(msr)THEN
340  DO i=1,elo_tm%NTIMES
341  READ(intelel,'(I7,300f8.4)') itemp,(rtemp1(j,i),j=1,ntidenode_gl)
342  END DO
343  END IF
344  CLOSE(intelel)
345 
346 
347 !
348 !---Map to Local Domain----------------------------------------
349 !
350  IF(serial) THEN
353  ALLOCATE(i_tidenode_n(ntidenode))
355  ALLOCATE(eltts(ntidenode,elo_tm%NTIMES)); eltts = zero
356  eltts = rtemp1
357  ENDIF
358 
359 
360  DEALLOCATE(rtemp1)
361 
362 
363  ELSE ! if statement end for ntidenode_GL > 0
364  close(intnode)
365  END IF
366 
367  CALL read_meanflow
368  CALL set_bndry_meanflow
369 
370 !
371 !--- calculate mapping function
372 !
373  ALLOCATE(i_obc_cell(0:nt),i_obc_node(0:mt),i_obc_cell2(0:nt),i_obc_node2(0:mt))
374 
375  i_obc_node = 0
376 ! I_OBC_NODE is the mapping from local domain node index j(MT) to tidal open
377 ! bndy node index i. If not a tidal open bndy node, I_OBC_NODE(j)=0
378  if(ntidenode > 0)then
379  do j = 1, mt
380  do i = 1, ntidenode
381  i1 = i_tidenode_n(i)
382  if (i1 == j) then
383  i_obc_node(j) = i
384  endif
385  enddo
386  enddo
387  endif
388 
389  i_obc_node2 = 0
390 ! I_OBC_NODE2 is the mapping from local domain node index j(MT) to open
391 ! bndy node index i. If not a open bndy node, I_OBC_NODE2(j)=0
392 !-------------------------------Jianzhong----------------------------!
393 ! if(iobcn > 0)then
394 ! do j = 1, MT
395 ! do i = 1, iobcn
396 ! I1 = I_OBC_N(i)
397 ! if (I1 == j) then
398 ! I_OBC_NODE2(J) = i
399 ! endif
400 ! enddo
401 ! enddo
402 ! endif
403  if(ibcn(2) > 0)then
404  do j = 1, mt
405  do i = 1, ibcn(2)
406  jn=obc_lst(2,i)
407  i1 = i_obc_n(jn)
408  if (i1 == j) then
409  i_obc_node2(j) = i
410  endif
411  enddo
412  enddo
413  endif
414 !--------------------------------------------------------------------!
415 
416  i_obc_cell = 0
417 ! I_OBC_CELL is the mapping from local domain cell index j(NT) to tidal open
418 ! bndy cell index i. if not a tidal open bndy cell, I_OBC_CELL(j)=0
419  if(ntidecell > 0)then
420  do j = 1, nt
421  do i = 1, ntidecell
422  i1 = i_tidecell_n(i)
423  if (i1 == j) then
424  i_obc_cell(j) = i
425  endif
426  enddo
427  enddo
428  endif
429 
430  i_obc_cell2 = 0
431 ! I_OBC_CELL2 is the mapping from local domain cell index j(NT) to mean flow open
432 ! bndy cell index i. if not a mean flow open bndy cell, I_OBC_CELL2(j)=0
433  if(nmfcell > 0)then
434  do j = 1, nt
435  do i = 1, nmfcell
436  i1 = i_mfcell_n(i)
437  if (i1 == j) then
438  i_obc_cell2(j) = i
439  endif
440  enddo
441  enddo
442  endif
443 
444 
445  IF(nmfcell > 0)THEN
446  ALLOCATE (dltn(nmfcell))
447  DO i = 1,nmfcell
448  ii = i_mfcell_n(i)
449  IF(nbe(ii,1) == 0 .and. isonb(nv(ii,1)) /= 2) THEN
450  deltx = vx(nv(ii,2))-vx(nv(ii,3))
451  delty = vy(nv(ii,2))-vy(nv(ii,3))
452  aa1 = -delty
453  bb1 = deltx
454  cc1 = -aa1*vx(nv(ii,2))-bb1*vy(nv(ii,2))
455  ELSE IF(nbe(ii,2) == 0 .and. isonb(nv(ii,2)) /= 2) THEN
456  deltx = vx(nv(ii,3))-vx(nv(ii,1))
457  delty = vy(nv(ii,3))-vy(nv(ii,1))
458  aa1 = -delty
459  bb1 = deltx
460  cc1 = -aa1*vx(nv(ii,3))-bb1*vy(nv(ii,3))
461  ELSE IF(nbe(ii,3) == 0 .and. isonb(nv(ii,3)) /= 2) THEN
462  deltx = vx(nv(ii,1))-vx(nv(ii,2))
463  delty = vy(nv(ii,1))-vy(nv(ii,2))
464  aa1 = -delty
465  bb1 = deltx
466  cc1 = -aa1*vx(nv(ii,1))-bb1*vy(nv(ii,1))
467  ELSE
468  print*,'something is wrong here 1'
469  CALL pstop
470  END IF
471 
472  aa2 = bb1
473  bb2 = -aa1
474  cc2 = -aa2*xc(ii)-bb2*yc(ii)
475 
476  xtmp1 = -(cc1*bb2-cc2*bb1)/(aa1*bb2-aa2*bb1)
477  ytmp1 = -(cc1*aa2-cc2*aa1)/(bb1*aa2-bb2*aa1)
478 
479  dltn(i) = sqrt((xc(ii)-xtmp1)**2+(yc(ii)-ytmp1)**2)
480  END DO
481  END IF
482 
483  RETURN
484  END SUBROUTINE find_obside
485 
486 
487 !==============================================================================|
488 ! INTERPOLATION of VARIOUS TIME SERIES |
489 !==============================================================================|
490 
491  SUBROUTINE bcond_tide_2d
493  USE bcs
494  USE mod_obcs
495  INTEGER L1,L2,IERR
496  REAL(SP) :: FACT,UFACT,TIME1
497  REAL(SP) :: TIMERK
498 
499  timerk = seconds(rktime-starttime)
500 ! TIME1 = TIMERK * 86400.0_SP
501  time1 = timerk
502 
503  CALL bracket(elo_tm,time1,l1,l2,fact,ufact,ierr)
504  IF(ntidecell > 0)THEN
505  uatf(1:ntidecell) = ufact*uatts(1:ntidecell,l1) + fact*uatts(1:ntidecell,l2)
506  vatf(1:ntidecell) = ufact*vatts(1:ntidecell,l1) + fact*vatts(1:ntidecell,l2)
507  uatf(1:ntidecell) = uatf(1:ntidecell) * ramp
508  vatf(1:ntidecell) = vatf(1:ntidecell) * ramp
509  END IF
510 
511  IF(ntidenode > 0)THEN
512  eltf(1:ntidenode) = ufact*eltts(1:ntidenode,l1) + fact*eltts(1:ntidenode,l2)
513  eltf(1:ntidenode) = eltf(1:ntidenode) * ramp
514  END IF
515 
516  RETURN
517  END SUBROUTINE bcond_tide_2d
518 
519 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
520 
521  SUBROUTINE bcond_tide_3d
522  USE bcs
523  USE mod_obcs
524  INTEGER L1,L2,IERR
525  REAL(SP) :: FACT,UFACT,TIME1
526 
527 ! TIME1 = TIME * 86400.0_SP - dti ! pay attention to this TIME (different between 2D and 3D)
528  time1=seconds(inttime-starttime)-dti
529  IF(ntidecell > 0)THEN
530  CALL bracket(elo_tm,time1,l1,l2,fact,ufact,ierr)
531  ut(1:ntidecell,:) = ufact*utts(1:ntidecell,:,l1) + fact*utts(1:ntidecell,:,l2)
532  vt(1:ntidecell,:) = ufact*vtts(1:ntidecell,:,l1) + fact*vtts(1:ntidecell,:,l2)
533  ut(1:ntidecell,:) = ut(1:ntidecell,:) * ramp
534  vt(1:ntidecell,:) = vt(1:ntidecell,:) * ramp
535  END IF
536 
537  RETURN
538  END SUBROUTINE bcond_tide_3d
539 
540 
541 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
542 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
543 
544 !==========================================================================|
545  SUBROUTINE bcond_ng_2d
547 !--------------------------------------------------------------------------|
548 ! NON-Gradient Open Boundary Condition (2-D) |
549 !--------------------------------------------------------------------------|
550 
551  USE all_vars
552  IMPLICIT NONE
553 
554  INTEGER :: I,J,J1
555 
556  IF(ntidecell > 0)THEN
557  DO i = 1, ntidecell
558  j = i_tidecell_n(i)
559  uap(i) = ua(j) - uat(i)
560  vap(i) = va(j) - vat(i)
561  END DO
562  END IF
563 
564  IF(nmfcell > 0)THEN
565  DO i = 1, nmfcell
566  j = i_mfcell_n(i)
567  j1= i_obc_cell(j)
568  uant(i) = uat(j1)
569  vant(i) = vat(j1)
570  uan(i) = ua(j )
571  van(i) = va(j )
572  uanp(i) = uap(j1)
573  vanp(i) = vap(j1)
574  END DO
575  END IF
576 
577  RETURN
578  END SUBROUTINE bcond_ng_2d
579 !==========================================================================|
580 
581 !==========================================================================|
582  SUBROUTINE bcond_ng_3d
584 !--------------------------------------------------------------------------|
585 ! NON-Gradient Open Boundary Condition (3-D) |
586 !--------------------------------------------------------------------------|
587 
588  USE all_vars
589  IMPLICIT NONE
590 
591  INTEGER :: I,J,J1,K
592 
593  IF(nmfcell > 0)THEN
594  DO i = 1, nmfcell
595  j = i_mfcell_n(i)
596  j1= i_obc_cell(j)
597  DO k = 1,kbm1
598  unt(i,k) = ut(j1,k)
599  vnt(i,k) = vt(j1,k)
600  un(i,k) = u(j ,k)
601  vn(i,k) = v(j ,k)
602  END DO
603  END DO
604  END IF
605 
606  RETURN
607  END SUBROUTINE bcond_ng_3d
608 !==========================================================================|
609 
610 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
611 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
612 
613 !==========================================================================|
614  SUBROUTINE bcond_bki_2d(KTT)
616 !--------------------------------------------------------------------------|
617 ! BLUMBERG AND KANTHA IMPLICIT OPEN BOUNDARY CONDITION (BKI) |
618 !--------------------------------------------------------------------------|
619 
620  USE all_vars
621  IMPLICIT NONE
622 
623  INTEGER, INTENT(IN) :: KTT
624  INTEGER :: I,II,J1,J
625  REAL(SP):: CC,CP
626  REAL(SP):: coef
627 
628 
629  coef = 10800.00_sp
630 
631  IF(ntidecell > 0)THEN
632  DO i = 1, ntidecell
633  j = i_tidecell_n(i)
634  uapf(i) = uaf(j) - uatf(i)
635  vapf(i) = vaf(j) - vatf(i)
636  uap(i) = uapf(i)
637  vap(i) = vapf(i)
638  END DO
639  END IF
640 
641  IF(nmfcell > 0)THEN
642  DO i = 1,nmfcell
643  ii = i_mfcell_n(i)
644  j1= i_obc_cell(ii)
645 
646  cc = sqrt(grav_e(ii)*d1(ii))*alpha_rk(ktt)*dte/dltn(i)
647  cp = cc + 1.0_sp
648  uantf(i) = (cc*uatf(j1) + uarknt(i)*(1.0_sp-alpha_rk(ktt)*dte/coef))/cp
649  vantf(i) = (cc*vatf(j1) + varknt(i)*(1.0_sp-alpha_rk(ktt)*dte/coef))/cp
650 ! UANF (I) = (CC*UAF (II) + UARKN(I) *(1.0_SP-ALPHA_RK(KTT)*DTE/coef))/CP
651 ! VANF (I) = (CC*VAF (II) + VARKN(I) *(1.0_SP-ALPHA_RK(KTT)*DTE/coef))/CP
652  uanf(i) = 0.0_sp
653  vanf(i) = 0.0_sp
654 
655  uanpf(i) = uanf(i) - uantf(i)
656  vanpf(i) = vanf(i) - vantf(i)
657  uant(i) = uantf(i)
658  vant(i) = vantf(i)
659  uan(i) = uanf(i)
660  van(i) = vanf(i)
661  uanp(i) = uanpf(i)
662  vanp(i) = vanpf(i)
663  END DO
664  END IF
665 
666  RETURN
667  END SUBROUTINE bcond_bki_2d
668 !==========================================================================|
669 
670 !==========================================================================|
671  SUBROUTINE bcond_bki_3d(KTT)
673 !--------------------------------------------------------------------------|
674 ! BLUMBERG AND KANTHA IMPLICIT OPEN BOUNDARY CONDITION (BKI) |
675 !--------------------------------------------------------------------------|
676 
677  USE all_vars
678  IMPLICIT NONE
679 
680  INTEGER, INTENT(IN) :: KTT
681  INTEGER :: I,II,K,J1
682  REAL(SP):: CC,CP
683  REAL(SP):: coef
684 
685  coef = 10800.00_sp
686 
687  IF(nmfcell > 0)THEN
688  DO i = 1,nmfcell
689  ii = i_mfcell_n(i)
690  j1= i_obc_cell(ii)
691 
692  cc = sqrt(grav_e(ii)*d1(ii))*dti/dltn(i)
693  cp = cc + 1.0_sp
694  DO k=1,kbm1
695  unt(i,k) = (cc*ut(j1,k) + untb(i,k)*(1.0_sp-dti/coef))/cp
696  vnt(i,k) = (cc*vt(j1,k) + vntb(i,k)*(1.0_sp-dti/coef))/cp
697 ! UN (I,K) = (CC*U (II,K) + UNB (I,K)*(1.0_SP-DTI/coef))/CP
698 ! VN (I,K) = (CC*V (II,K) + VNB (I,K)*(1.0_SP-DTI/coef))/CP
699  un(i,k) = 0.0_sp
700  vn(i,k) = 0.0_sp
701  IF (ktt == 2) THEN
702  untb(i,k) = unt(i,k)
703  vntb(i,k) = vnt(i,k)
704  unb(i,k) = un(i,k)
705  vnb(i,k) = vn(i,k)
706  END IF
707  END DO
708  END DO
709  END IF
710 
711  RETURN
712  END SUBROUTINE bcond_bki_3d
713 !==========================================================================|
714 
715 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
716 !--------------------------------------------------------------------------|
717 ! TEST WHETHER THE TWO BOUNDARY CONDITION WORK |
718 !--------------------------------------------------------------------------|
719  SUBROUTINE test_cell(INDEX,MSG)
721 
722  IMPLICIT NONE
723 
724  CHARACTER(LEN=*) :: MSG
725  INTEGER :: INDEX
726  INTEGER,PARAMETER :: CELL1=59 ! IN TAIWAN STRAIT
727  INTEGER,PARAMETER :: CELL2=229 ! IN SECOND OPEN BOUNDARY
728 
729 
730 
731 ! IF(CELL1==EGID(INDEX))THEN
732 ! WRITE(IPT_P,*)'IN '//MSG//' CELL OF TAIWAN STARIT IN',MYID
733 ! write(ipt_p,*) "INDEX:",index,"egid(index)",egid(index),"; msg="//TRIM(MSG)//"; CELLS:",CELL1,CELL2
734 ! END IF
735 ! IF(CELL2==EGID(INDEX))THEN
736 ! WRITE(IPT_P,*)'IN '//MSG//' CELL OF PACIFIC BOUNDARY IN',MYID
737 ! write(ipt_p,*) "INDEX:",index,"egid(index)",egid(index),"; msg="//TRIM(MSG)//"; CELLS:",CELL1,CELL2
738 ! END IF
739 
740  RETURN
741  END SUBROUTINE test_cell
742 !==========================================================================|
743 
744  SUBROUTINE test_node(INDEX,MSG)
746  IMPLICIT NONE
747 
748  CHARACTER(LEN=*) :: MSG
749  INTEGER :: INDEX
750  INTEGER,PARAMETER :: NODE1=30 ! IN TAIWAN STRAIT
751  INTEGER,PARAMETER :: NODE2=115 ! IN SECOND OPEN BOUNDARY
752 
753 
754 ! IF(NODE1==NGID(INDEX))THEN
755 ! WRITE(IPT_P,*)'IN '//MSG//' NODE OF TAIWAN STARIT IN',MYID
756 ! write(ipt_p,*) "INDEX:",index,"ngid(index)",NGID(index),"; msg="//TRIM(MSG)//"; NODES:",NODE1,NODE2
757 ! END IF
758 ! IF(NODE2==NGID(INDEX))THEN
759 ! WRITE(IPT_P,*)'IN '//MSG//' NODE OF PACIFIC BOUNDARY IN',MYID
760 ! write(ipt_p,*) "INDEX:",index,"ngid(index)",NGID(index),"; msg="//TRIM(MSG)//"; NODES:",NODE1,NODE2
761 ! END IF
762 
763  RETURN
764  END SUBROUTINE test_node
765 !==========================================================================|
766 
767 
768 !%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%|
769 
770 END MODULE mod_obcs2
real(sp), dimension(:), allocatable, target va
Definition: mod_main.f90:1104
real(sp), dimension(:,:), allocatable untb
Definition: mod_obcs2.f90:73
real(sp), dimension(:), allocatable eltdt
Definition: mod_obcs2.f90:58
type(bc) elo_tm
Definition: mod_main.f90:1767
subroutine bcond_bki_2d(KTT)
Definition: mod_obcs2.f90:615
real(sp), dimension(:), allocatable, target d1
Definition: mod_main.f90:1116
real(sp), dimension(:), allocatable uapf
Definition: mod_obcs2.f90:69
integer, dimension(:), allocatable i_obc_cell
Definition: mod_obcs2.f90:55
real(sp), dimension(:), allocatable vanp
Definition: mod_obcs2.f90:64
subroutine test_cell(INDEX, MSG)
Definition: mod_obcs2.f90:720
integer ntidenode
Definition: mod_obcs2.f90:76
real(sp), dimension(:,:), allocatable, target v
Definition: mod_main.f90:1269
real(sp), dimension(:,:), allocatable vntb
Definition: mod_obcs2.f90:73
real(sp), dimension(:,:,:), allocatable utts
Definition: mod_obcs2.f90:56
integer, dimension(:), allocatable i_mfcell_n
real(sp), dimension(:), allocatable uan
Definition: mod_obcs2.f90:63
subroutine bcond_ng_2d
Definition: mod_obcs2.f90:546
real(sp), dimension(:), allocatable van
Definition: mod_obcs2.f90:63
integer, dimension(:), allocatable i_tidenode_n
Definition: mod_obcs2.f90:77
real(sp), dimension(:), allocatable vapf
Definition: mod_obcs2.f90:69
integer intelel
real(sp), dimension(:), allocatable, target yc
Definition: mod_main.f90:1004
real(sp), dimension(:), allocatable elp
Definition: mod_obcs2.f90:59
integer, dimension(:), allocatable i_obc_node
Definition: mod_obcs2.f90:55
real(sp), dimension(:,:), allocatable eltts
Definition: mod_obcs2.f90:56
integer, dimension(:,:), allocatable obc_lst
Definition: mod_obcs.f90:84
real(sp), dimension(:), allocatable varkn
Definition: mod_obcs2.f90:72
integer ntidenode_i
Definition: mod_obcs2.f90:76
subroutine test_node(INDEX, MSG)
Definition: mod_obcs2.f90:745
type(time) starttime
Definition: mod_main.f90:833
real(sp), dimension(:), allocatable uanp
Definition: mod_obcs2.f90:64
real(sp), dimension(:,:), allocatable un
Definition: mod_obcs2.f90:67
real(sp), dimension(:), allocatable uanf
Definition: mod_obcs2.f90:70
real(sp), dimension(:), allocatable elrkt
Definition: mod_obcs2.f90:58
real(sp), dimension(:,:), allocatable ut
Definition: mod_obcs2.f90:66
real(sp), dimension(:,:), allocatable uatts
Definition: mod_obcs2.f90:56
real(sp), dimension(:), allocatable dltn
Definition: mod_obcs2.f90:80
real(sp), dimension(:,:), allocatable, target u
Definition: mod_main.f90:1268
real(sp), dimension(:), allocatable uantf
Definition: mod_obcs2.f90:70
integer, dimension(:,:), allocatable, target iec
Definition: mod_main.f90:1028
real(sp), dimension(:), allocatable elt
Definition: mod_obcs2.f90:58
real(sp), dimension(:), allocatable uatf
Definition: mod_obcs2.f90:60
real(sp), dimension(:), allocatable uarknt
Definition: mod_obcs2.f90:72
subroutine find_obside
Definition: mod_obcs2.f90:138
integer, dimension(:,:), allocatable, target niec
Definition: mod_main.f90:1032
subroutine bcond_ng_3d
Definition: mod_obcs2.f90:583
real(sp), dimension(:,:), allocatable vn
Definition: mod_obcs2.f90:67
integer, dimension(:), allocatable cobedge_lst
Definition: mod_obcs2.f90:54
integer ntidecell
Definition: mod_obcs2.f90:75
real(sp), dimension(:), allocatable uap
Definition: mod_obcs2.f90:61
real(sp), dimension(:,:), allocatable vnt
Definition: mod_obcs2.f90:67
real(sp), dimension(:), allocatable, target vx
Definition: mod_main.f90:1001
real(sp), dimension(:), allocatable, target grav_e
Definition: mod_main.f90:1013
integer ntidenode_gl
Definition: mod_obcs2.f90:76
real(sp), dimension(:), allocatable, target vaf
Definition: mod_main.f90:1106
integer, parameter sp
Definition: mod_prec.f90:48
real(sp), dimension(:), allocatable vanpf
Definition: mod_obcs2.f90:71
real(sp), dimension(:), allocatable, target vy
Definition: mod_main.f90:1002
integer, dimension(:), allocatable i_tidecell_gl
Definition: mod_obcs2.f90:78
integer, dimension(:), allocatable, target ntve
Definition: mod_main.f90:1022
real(sp), dimension(:), allocatable uarkn
Definition: mod_obcs2.f90:72
real(sp), dimension(:), allocatable uanpf
Definition: mod_obcs2.f90:71
integer, dimension(:,:), allocatable, target nbe
Definition: mod_main.f90:1020
integer, dimension(:), allocatable i_obc_n
Definition: mod_main.f90:1779
real(sp), dimension(:), allocatable vantf
Definition: mod_obcs2.f90:70
integer intnode
subroutine bcond_bki_3d(KTT)
Definition: mod_obcs2.f90:672
integer nobcv
Definition: mod_obcs2.f90:53
subroutine alloc_obc2_data
Definition: mod_obcs2.f90:86
integer, dimension(:,:), allocatable, target nv
Definition: mod_main.f90:1018
integer, dimension(5) ibcn
Definition: mod_obcs.f90:82
real(sp), dimension(:), allocatable vanf
Definition: mod_obcs2.f90:70
integer nmfcell
real(sp), dimension(:,:), allocatable vt
Definition: mod_obcs2.f90:66
real(sp), dimension(:), allocatable uant
Definition: mod_obcs2.f90:63
subroutine bcond_tide_3d
Definition: mod_obcs2.f90:522
real(sp), dimension(:,:), allocatable unt
Definition: mod_obcs2.f90:67
real(sp), dimension(:), allocatable elrkp
Definition: mod_obcs2.f90:59
real(sp), dimension(:), allocatable, target ua
Definition: mod_main.f90:1103
real(sp), dimension(:), allocatable vant
Definition: mod_obcs2.f90:63
real(sp), dimension(:,:), allocatable vatts
Definition: mod_obcs2.f90:56
subroutine read_meanflow
real(sp), dimension(:,:), allocatable unb
Definition: mod_obcs2.f90:73
real(sp), dimension(:), allocatable uat
Definition: mod_obcs2.f90:60
real(sp), dimension(:), allocatable elpf
Definition: mod_obcs2.f90:59
integer, dimension(:,:), allocatable, target nbve
Definition: mod_main.f90:1034
integer, dimension(:), allocatable i_obc_cell2
Definition: mod_obcs2.f90:55
integer, dimension(:), allocatable i_obc_node2
Definition: mod_obcs2.f90:55
real(sp), dimension(:), allocatable, target xc
Definition: mod_main.f90:1003
real(sp), dimension(:), allocatable, target uaf
Definition: mod_main.f90:1105
real(sp), dimension(:), allocatable eltf
Definition: mod_obcs2.f90:58
integer, dimension(:), allocatable nobedge_lst
Definition: mod_obcs2.f90:54
real(sp) ramp
Definition: mod_main.f90:845
integer ntidecell_gl
Definition: mod_obcs2.f90:75
integer ntidecell_i
Definition: mod_obcs2.f90:75
real(sp), dimension(:,:), allocatable vnb
Definition: mod_obcs2.f90:73
subroutine bracket(TMAP, STIME, L1, L2, FACT, BACT, IERR)
real(sp), dimension(:), allocatable varknt
Definition: mod_obcs2.f90:72
real(sp), dimension(:), allocatable vat
Definition: mod_obcs2.f90:60
integer nobe
Definition: mod_obcs2.f90:53
integer intcell
real(sp), dimension(:,:,:), allocatable vtts
Definition: mod_obcs2.f90:56
subroutine set_bndry_meanflow
integer, dimension(:), allocatable, target isonb
Definition: mod_main.f90:1024
integer ios
Definition: mod_obcs2.f90:81
integer, dimension(:), allocatable i_tidenode_gl
Definition: mod_obcs2.f90:77
subroutine bcond_tide_2d
Definition: mod_obcs2.f90:492
real(sp), dimension(:), allocatable vap
Definition: mod_obcs2.f90:61
integer, dimension(:), allocatable i_tidecell_n
Definition: mod_obcs2.f90:78
real(sp), dimension(:), allocatable vatf
Definition: mod_obcs2.f90:60
type(time) rktime
Definition: mod_main.f90:829