Back to home page

MITgcm

 
 

    


File indexing completed on 2024-03-19 05:10:11 UTC

view on githubraw file Latest commit 720a211d on 2024-03-18 17:10:23 UTC
7bfe6112e8 Jean*0001 #include "CTRL_OPTIONS.h"
7109a141b2 Patr*0002 #ifdef ALLOW_OBCS
                0003 # include "OBCS_OPTIONS.h"
                0004 #endif
                0005 
                0006       subroutine ctrl_getobcsn(
9f5240b52a Jean*0007      I                          myTime, myIter, myThid )
7109a141b2 Patr*0008 
                0009 c     ==================================================================
                0010 c     SUBROUTINE ctrl_getobcsn
                0011 c     ==================================================================
                0012 c
                0013 c     o Get northern obc of the control vector and add it
                0014 c       to dyn. fields
                0015 c
                0016 c     started: heimbach@mit.edu, 29-Aug-2001
                0017 c     modified: gebbie@mit.edu, 18-Mar-2003
9f5240b52a Jean*0018 c
7109a141b2 Patr*0019 c     ==================================================================
                0020 c     SUBROUTINE ctrl_getobcsn
                0021 c     ==================================================================
                0022 
                0023       implicit none
                0024 
                0025 c     == global variables ==
46c1d12550 Jean*0026 #ifdef ALLOW_OBCSN_CONTROL
7109a141b2 Patr*0027 #include "EEPARAMS.h"
                0028 #include "SIZE.h"
                0029 #include "PARAMS.h"
                0030 #include "GRID.h"
46c1d12550 Jean*0031 c#include "OBCS_PARAMS.h"
                0032 #include "OBCS_GRID.h"
                0033 #include "OBCS_FIELDS.h"
c04085ad02 Patr*0034 #include "CTRL_SIZE.h"
4d72283393 Mart*0035 #include "CTRL.h"
edcd27be69 Mart*0036 #include "CTRL_DUMMY.h"
e612621177 Gael*0037 #include "CTRL_OBCS.h"
65754df434 Mart*0038 #include "OPTIMCYCLE.h"
46c1d12550 Jean*0039 #endif /* ALLOW_OBCSN_CONTROL */
7109a141b2 Patr*0040 
                0041 c     == routine arguments ==
9f5240b52a Jean*0042       _RL     myTime
                0043       integer myIter
                0044       integer myThid
7109a141b2 Patr*0045 
46c1d12550 Jean*0046 #ifdef ALLOW_OBCSN_CONTROL
9f5240b52a Jean*0047 c     == external functions ==
                0048       integer  ilnblnk
                0049       external ilnblnk
7109a141b2 Patr*0050 
9f5240b52a Jean*0051 c     == local variables ==
7109a141b2 Patr*0052       integer bi,bj
                0053       integer i,j,k
                0054       integer itlo,ithi
                0055       integer jtlo,jthi
                0056       integer imin,imax
                0057       integer ilobcsn
                0058       integer iobcs
                0059       integer jp1
                0060       _RL     obcsnfac
                0061       logical obcsnfirst
                0062       logical obcsnchanged
                0063       integer obcsncount0
                0064       integer obcsncount1
9f5240b52a Jean*0065 cgg      _RL maskxz   (1-OLx:sNx+OLx,Nr,nSx,nSy)
                0066       _RL tmpfldxz (1-OLx:sNx+OLx,Nr,nSx,nSy)
7109a141b2 Patr*0067       logical doglobalread
                0068       logical ladinit
de57a2ec4b Mart*0069       character*(MAX_LEN_FNAM) fnameobcsn
b6199164d9 Matt*0070 #ifdef ALLOW_OBCS_CONTROL_MODES
                0071       integer nk,nz
9f5240b52a Jean*0072       _RL     tmpz (Nr,nSx,nSy)
b6199164d9 Matt*0073       _RL     stmp
                0074 #endif
f9d7cbfb72 Ou W*0075       integer ilDir
7109a141b2 Patr*0076 
                0077 c     == end of interface ==
                0078 
9f5240b52a Jean*0079       jtlo = myByLo(myThid)
                0080       jthi = myByHi(myThid)
                0081       itlo = myBxLo(myThid)
                0082       ithi = myBxHi(myThid)
7109a141b2 Patr*0083 cgg      imin = 1-olx
                0084 cgg      imax = snx+olx
                0085       imin = 1
9f5240b52a Jean*0086       imax = sNx
7109a141b2 Patr*0087       jp1  = 0
                0088 
                0089 c--   Now, read the control vector.
                0090       doglobalread = .false.
                0091       ladinit      = .false.
                0092 
f9d7cbfb72 Ou W*0093 c     Find ctrlDir (w/o trailing blanks) length
                0094       ilDir = ilnblnk(ctrlDir)
                0095 
7109a141b2 Patr*0096       if (optimcycle .ge. 0) then
                0097         ilobcsn=ilnblnk( xx_obcsn_file )
de57a2ec4b Mart*0098         write(fnameobcsn,'(2a,i10.10)')
f9d7cbfb72 Ou W*0099      &       ctrlDir(1:ilDir)//xx_obcsn_file(1:ilobcsn), '.', optimcycle
7109a141b2 Patr*0100       endif
                0101 
                0102 c--   Get the counters, flags, and the interpolation factor.
                0103       call ctrl_get_gen_rec(
                0104      I                   xx_obcsnstartdate, xx_obcsnperiod,
                0105      O                   obcsnfac, obcsnfirst, obcsnchanged,
                0106      O                   obcsncount0,obcsncount1,
9f5240b52a Jean*0107      I                   myTime, myIter, myThid )
7109a141b2 Patr*0108 
                0109       do iobcs = 1,nobcs
9f5240b52a Jean*0110 
859d62c7fb Mart*0111        if ( obcsnfirst ) then
720a211d38 Ou W*0112 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0113         call active_read_xz( fnameobcsn, tmpfldxz,
                0114      &                       (obcsncount0-1)*nobcs+iobcs,
                0115      &                       doglobalread, ladinit, optimcycle,
4d72283393 Mart*0116      &                       myThid, xx_obcsn_dummy )
720a211d38 Ou W*0117 #else
                0118         CALL READ_REC_XZ_RL( fnameobcsn, ctrlprec, Nr, tmpfldxz,
                0119      &                       (obcsncount0-1)*nobcs+iobcs, 1, myThid )
                0120 #endif
859d62c7fb Mart*0121 
                0122         do bj = jtlo,jthi
                0123          do bi = itlo,ithi
b6199164d9 Matt*0124 #ifdef ALLOW_OBCS_CONTROL_MODES
                0125           if (iobcs .gt. 2) then
                0126            do i = imin,imax
                0127             j = OB_Jn(i,bi,bj)
9fd29e65a3 Jean*0128             IF ( j.EQ.OB_indexNone ) j = 1
b6199164d9 Matt*0129 cih    Determine number of open vertical layers.
                0130             nz = 0
                0131             do k = 1,Nr
                0132              if (iobcs .eq. 3) then
                0133               nz = nz + maskS(i,j+jp1,k,bi,bj)
                0134              else
                0135               nz = nz + maskW(i,j,k,bi,bj)
                0136              endif
                0137             end do
                0138 cih    Compute absolute velocities from the barotropic-baroclinic modes.
                0139             do k = 1,Nr
                0140              if (k.le.nz) then
                0141               stmp = 0.
                0142               do nk = 1,nz
                0143                stmp = stmp +
                0144      &         modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
                0145               end do
                0146               tmpz(k,bi,bj) = stmp
                0147              else
                0148               tmpz(k,bi,bj) = 0.
                0149              end if
                0150             end do
                0151             do k = 1,Nr
                0152              if (iobcs .eq. 3) then
                0153               tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
                0154      &         *recip_hFacS(i,j+jp1,k,bi,bj)
                0155              else
                0156               tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
                0157      &                           *recip_hFacW(i,j,k,bi,bj)
                0158              endif
                0159             end do
                0160            enddo
                0161           endif
                0162 #endif
9f5240b52a Jean*0163           do k = 1,Nr
859d62c7fb Mart*0164            do i = imin,imax
                0165             xx_obcsn1(i,k,bi,bj,iobcs)  = tmpfldxz (i,k,bi,bj)
                0166 cgg   &                                    *   maskxz (i,k,bi,bj)
7c7f2df94e Mart*0167            enddo
7109a141b2 Patr*0168           enddo
859d62c7fb Mart*0169          enddo
                0170         enddo
                0171        endif
7109a141b2 Patr*0172 
859d62c7fb Mart*0173        if ( (obcsnfirst) .or. (obcsnchanged)) then
46c1d12550 Jean*0174 
859d62c7fb Mart*0175         do bj = jtlo,jthi
                0176          do bi = itlo,ithi
9f5240b52a Jean*0177           do k = 1,Nr
859d62c7fb Mart*0178            do i = imin,imax
                0179             xx_obcsn0(i,k,bi,bj,iobcs) = xx_obcsn1(i,k,bi,bj,iobcs)
                0180             tmpfldxz (i,k,bi,bj)       = 0. _d 0
                0181            enddo
7109a141b2 Patr*0182           enddo
859d62c7fb Mart*0183          enddo
                0184         enddo
46c1d12550 Jean*0185 
720a211d38 Ou W*0186 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0187         call active_read_xz( fnameobcsn, tmpfldxz,
                0188      &                       (obcsncount1-1)*nobcs+iobcs,
                0189      &                       doglobalread, ladinit, optimcycle,
4d72283393 Mart*0190      &                       myThid, xx_obcsn_dummy )
720a211d38 Ou W*0191 #else
                0192         CALL READ_REC_XZ_RL( fnameobcsn, ctrlprec, Nr, tmpfldxz,
                0193      &                       (obcsncount1-1)*nobcs+iobcs, 1, myThid )
                0194 #endif
7109a141b2 Patr*0195 
                0196         do bj = jtlo,jthi
7aa90384e1 Mart*0197          do bi = itlo,ithi
b6199164d9 Matt*0198 #ifdef ALLOW_OBCS_CONTROL_MODES
                0199           if (iobcs .gt. 2) then
                0200            do i = imin,imax
                0201             j = OB_Jn(i,bi,bj)
9fd29e65a3 Jean*0202             IF ( j.EQ.OB_indexNone ) j = 1
b6199164d9 Matt*0203 cih    Determine number of open vertical layers.
                0204             nz = 0
                0205             do k = 1,Nr
                0206              if (iobcs .eq. 3) then
                0207               nz = nz + maskS(i,j+jp1,k,bi,bj)
                0208              else
                0209               nz = nz + maskW(i,j,k,bi,bj)
                0210              endif
                0211             end do
                0212 cih    Compute absolute velocities from the barotropic-baroclinic modes.
                0213             do k = 1,Nr
                0214              if (k.le.nz) then
                0215               stmp = 0.
                0216               do nk = 1,nz
                0217                stmp = stmp +
                0218      &         modesv(k,nk,nz)*tmpfldxz(i,nk,bi,bj)
                0219               end do
                0220               tmpz(k,bi,bj) = stmp
                0221              else
                0222               tmpz(k,bi,bj) = 0.
                0223              end if
                0224             end do
                0225             do k = 1,Nr
                0226              if (iobcs .eq. 3) then
                0227               tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
                0228      &         *recip_hFacS(i,j+jp1,k,bi,bj)
                0229              else
                0230               tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
                0231      &                           *recip_hFacW(i,j,k,bi,bj)
                0232              endif
                0233             end do
                0234            enddo
                0235           endif
                0236 #endif
9f5240b52a Jean*0237           do k = 1,Nr
859d62c7fb Mart*0238            do i = imin,imax
                0239             xx_obcsn1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
                0240 cgg   &                                        *   maskxz (i,k,bi,bj)
7109a141b2 Patr*0241            enddo
7aa90384e1 Mart*0242           enddo
                0243          enddo
7109a141b2 Patr*0244         enddo
859d62c7fb Mart*0245        endif
                0246 
                0247 c--   Add control to model variable.
                0248        do bj = jtlo,jthi
                0249         do bi = itlo,ithi
                0250 c--   Calculate mask for tracer cells (0 => land, 1 => water).
9f5240b52a Jean*0251          do k = 1,Nr
                0252           do i = 1,sNx
                0253            j = OB_Jn(i,bi,bj)
9fd29e65a3 Jean*0254            IF ( j.EQ.OB_indexNone ) j = 1
859d62c7fb Mart*0255            if (iobcs .EQ. 1) then
                0256             OBNt(i,k,bi,bj) = OBNt (i,k,bi,bj)
                0257      &           + obcsnfac            *xx_obcsn0(i,k,bi,bj,iobcs)
                0258      &           + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
                0259             OBNt(i,k,bi,bj) = OBNt(i,k,bi,bj)
                0260      &           *maskS(i,j+jp1,k,bi,bj)
                0261            else if (iobcs .EQ. 2) then
                0262             OBNs(i,k,bi,bj) = OBNs (i,k,bi,bj)
                0263      &           + obcsnfac            *xx_obcsn0(i,k,bi,bj,iobcs)
                0264      &           + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
                0265             OBNs(i,k,bi,bj) = OBNs(i,k,bi,bj)
                0266      &           *maskS(i,j+jp1,k,bi,bj)
                0267            else if (iobcs .EQ. 4) then
                0268             OBNu(i,k,bi,bj) = OBNu (i,k,bi,bj)
                0269      &           + obcsnfac            *xx_obcsn0(i,k,bi,bj,iobcs)
                0270      &           + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
                0271             OBNu(i,k,bi,bj) = OBNu(i,k,bi,bj)
                0272      &           *maskW(i,j,k,bi,bj)
                0273            else if (iobcs .EQ. 3) then
                0274             OBNv(i,k,bi,bj) = OBNv (i,k,bi,bj)
                0275      &           + obcsnfac            *xx_obcsn0(i,k,bi,bj,iobcs)
                0276      &           + (1. _d 0 - obcsnfac)*xx_obcsn1(i,k,bi,bj,iobcs)
                0277             OBNv(i,k,bi,bj) = OBNv(i,k,bi,bj)
                0278      &           *maskS(i,j+jp1,k,bi,bj)
                0279            endif
                0280           enddo
                0281          enddo
                0282         enddo
                0283        enddo
46c1d12550 Jean*0284 
7109a141b2 Patr*0285 C--   End over iobcs loop
                0286       enddo
951926fb9b Jean*0287 
7109a141b2 Patr*0288 #endif /* ALLOW_OBCSN_CONTROL */
                0289 
46c1d12550 Jean*0290       return
7109a141b2 Patr*0291       end