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_getobcsw(
9f5240b52a Jean*0007      I                          myTime, myIter, myThid )
7109a141b2 Patr*0008 
                0009 c     ==================================================================
                0010 c     SUBROUTINE ctrl_getobcsw
                0011 c     ==================================================================
                0012 c
                0013 c     o Get western 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_getobcsw
                0021 c     ==================================================================
                0022 
                0023       implicit none
                0024 
                0025 c     == global variables ==
46c1d12550 Jean*0026 #ifdef ALLOW_OBCSW_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_OBCSW_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_OBCSW_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 jmin,jmax
                0057       integer ilobcsw
                0058       integer iobcs
                0059       integer ip1
                0060       _RL     obcswfac
                0061       logical obcswfirst
                0062       logical obcswchanged
                0063       integer obcswcount0
                0064       integer obcswcount1
9f5240b52a Jean*0065 cgg      _RL maskyz   (1-OLy:sNy+OLy,Nr,nSx,nSy)
                0066       _RL tmpfldyz (1-OLy:sNy+OLy,Nr,nSx,nSy)
7109a141b2 Patr*0067       logical doglobalread
                0068       logical ladinit
de57a2ec4b Mart*0069       character*(MAX_LEN_FNAM) fnameobcsw
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)
                0083       jmin = 1-OLy
                0084       jmax = sNy+OLy
7109a141b2 Patr*0085       ip1  = 1
                0086 
                0087 c--   Now, read the control vector.
                0088       doglobalread = .false.
                0089       ladinit      = .false.
                0090 
f9d7cbfb72 Ou W*0091 c     Find ctrlDir (w/o trailing blanks) length
                0092       ilDir = ilnblnk(ctrlDir)
                0093 
7109a141b2 Patr*0094       if (optimcycle .ge. 0) then
859d62c7fb Mart*0095        ilobcsw=ilnblnk( xx_obcsw_file )
de57a2ec4b Mart*0096        write(fnameobcsw,'(2a,i10.10)')
f9d7cbfb72 Ou W*0097      &      ctrlDir(1:ilDir)//xx_obcsw_file(1:ilobcsw), '.', optimcycle
7109a141b2 Patr*0098       endif
                0099 
                0100 c--   Get the counters, flags, and the interpolation factor.
                0101       call ctrl_get_gen_rec(
                0102      I                   xx_obcswstartdate, xx_obcswperiod,
                0103      O                   obcswfac, obcswfirst, obcswchanged,
                0104      O                   obcswcount0,obcswcount1,
9f5240b52a Jean*0105      I                   myTime, myIter, myThid )
7109a141b2 Patr*0106 
                0107       do iobcs = 1,nobcs
9f5240b52a Jean*0108 
859d62c7fb Mart*0109        if ( obcswfirst ) then
720a211d38 Ou W*0110 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0111         call active_read_yz( fnameobcsw, tmpfldyz,
                0112      &                       (obcswcount0-1)*nobcs+iobcs,
                0113      &                       doglobalread, ladinit, optimcycle,
4d72283393 Mart*0114      &                       myThid, xx_obcsw_dummy )
720a211d38 Ou W*0115 #else
                0116         CALL READ_REC_YZ_RL( fnameobcsw, ctrlprec, Nr, tmpfldyz,
                0117      &                       (obcswcount0-1)*nobcs+iobcs, 1, myThid )
                0118 #endif
859d62c7fb Mart*0119 
                0120         do bj = jtlo,jthi
                0121          do bi = itlo,ithi
b6199164d9 Matt*0122 #ifdef ALLOW_OBCS_CONTROL_MODES
                0123           if (iobcs .gt. 2) then
                0124            do j = jmin,jmax
                0125             i = OB_Iw(j,bi,bj)
9fd29e65a3 Jean*0126             IF ( i.EQ.OB_indexNone ) i = 1
b6199164d9 Matt*0127 cih    Determine number of open vertical layers.
                0128             nz = 0
                0129             do k = 1,Nr
                0130               if (iobcs .eq. 3) then
                0131                 nz = nz + maskW(i+ip1,j,k,bi,bj)
                0132               else
                0133                 nz = nz + maskS(i,j,k,bi,bj)
                0134               endif
                0135             end do
                0136 cih    Compute absolute velocities from the barotropic-baroclinic modes.
                0137             do k = 1,Nr
                0138              if (k.le.nz) then
                0139               stmp = 0.
                0140               do nk = 1,nz
                0141                stmp = stmp +
                0142      &         modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj)
                0143               end do
                0144                tmpz(k,bi,bj) = stmp
                0145              else
                0146               tmpz(k,bi,bj) = 0.
                0147              end if
                0148             enddo
                0149             do k = 1,Nr
                0150               if (iobcs .eq. 3) then
                0151                 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
                0152      &            *recip_hFacW(i+ip1,j,k,bi,bj)
                0153               else
                0154                 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
                0155      &            *recip_hFacS(i,j,k,bi,bj)
                0156               endif
                0157             end do
                0158            enddo
                0159           endif
                0160 #endif
9f5240b52a Jean*0161           do k = 1,Nr
859d62c7fb Mart*0162            do j = jmin,jmax
                0163             xx_obcsw1(j,k,bi,bj,iobcs)  = tmpfldyz (j,k,bi,bj)
                0164 cgg   &                                        *   maskyz (j,k,bi,bj)
7c7f2df94e Mart*0165            enddo
7109a141b2 Patr*0166           enddo
859d62c7fb Mart*0167          enddo
                0168         enddo
                0169        endif
                0170 
                0171        if ( (obcswfirst) .or. (obcswchanged)) then
7109a141b2 Patr*0172 
859d62c7fb Mart*0173         do bj = jtlo,jthi
                0174          do bi = itlo,ithi
9f5240b52a Jean*0175           do k = 1,Nr
859d62c7fb Mart*0176            do j = jmin,jmax
                0177             xx_obcsw0(j,k,bi,bj,iobcs) = xx_obcsw1(j,k,bi,bj,iobcs)
                0178             tmpfldyz (j,k,bi,bj)       = 0. _d 0
                0179            enddo
7109a141b2 Patr*0180           enddo
859d62c7fb Mart*0181          enddo
                0182         enddo
7109a141b2 Patr*0183 
720a211d38 Ou W*0184 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0185         call active_read_yz( fnameobcsw, tmpfldyz,
                0186      &                       (obcswcount1-1)*nobcs+iobcs,
                0187      &                       doglobalread, ladinit, optimcycle,
4d72283393 Mart*0188      &                       myThid, xx_obcsw_dummy )
720a211d38 Ou W*0189 #else
                0190         CALL READ_REC_YZ_RL( fnameobcsw, ctrlprec, Nr, tmpfldyz,
                0191      &                       (obcswcount1-1)*nobcs+iobcs, 1, myThid )
                0192 #endif
859d62c7fb Mart*0193 
                0194         do bj = jtlo,jthi
                0195          do bi = itlo,ithi
b6199164d9 Matt*0196 #ifdef ALLOW_OBCS_CONTROL_MODES
                0197           if (iobcs .gt. 2) then
                0198            do j = jmin,jmax
                0199             i = OB_Iw(j,bi,bj)
9fd29e65a3 Jean*0200             IF ( i.EQ.OB_indexNone ) i = 1
b6199164d9 Matt*0201 cih    Determine number of open vertical layers.
                0202             nz = 0
                0203             do k = 1,Nr
                0204               if (iobcs .eq. 3) then
                0205                 nz = nz + maskW(i+ip1,j,k,bi,bj)
                0206               else
                0207                 nz = nz + maskS(i,j,k,bi,bj)
                0208               endif
                0209             end do
                0210 cih    Compute absolute velocities from the barotropic-baroclinic modes.
                0211             do k = 1,Nr
                0212              if (k.le.nz) then
                0213               stmp = 0.
                0214               do nk = 1,nz
                0215                stmp = stmp +
                0216      &         modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj)
                0217               end do
                0218                tmpz(k,bi,bj) = stmp
                0219              else
                0220                tmpz(k,bi,bj) = 0.
                0221              end if
                0222             enddo
                0223             do k = 1,Nr
                0224               if (iobcs .eq. 3) then
                0225                 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
                0226      &            *recip_hFacW(i+ip1,j,k,bi,bj)
                0227               else
                0228                 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
                0229      &            *recip_hFacS(i,j,k,bi,bj)
                0230               endif
                0231             end do
                0232            enddo
                0233           endif
                0234 #endif
9f5240b52a Jean*0235           do k = 1,Nr
859d62c7fb Mart*0236            do j = jmin,jmax
                0237             xx_obcsw1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
                0238 cgg   &                                        *   maskyz (j,k,bi,bj)
7109a141b2 Patr*0239            enddo
7aa90384e1 Mart*0240           enddo
                0241          enddo
7109a141b2 Patr*0242         enddo
859d62c7fb Mart*0243        endif
7109a141b2 Patr*0244 
859d62c7fb Mart*0245 c--   Add control to model variable.
                0246        do bj = jtlo, jthi
                0247         do bi = itlo, ithi
                0248 c--   Calculate mask for tracer cells (0 => land, 1 => water).
9f5240b52a Jean*0249          do k = 1,Nr
                0250           do j = 1,sNy
859d62c7fb Mart*0251            i = OB_Iw(j,bi,bj)
9fd29e65a3 Jean*0252            IF ( i.EQ.OB_indexNone ) i = 1
859d62c7fb Mart*0253            if (iobcs .EQ. 1) then
                0254             OBWt(j,k,bi,bj) = OBWt (j,k,bi,bj)
                0255      &           + obcswfac            *xx_obcsw0(j,k,bi,bj,iobcs)
                0256      &           + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs)
                0257             OBWt(j,k,bi,bj) = OBWt(j,k,bi,bj)
                0258      &           *maskW(i+ip1,j,k,bi,bj)
                0259            else if (iobcs .EQ. 2) then
                0260             OBWs(j,k,bi,bj) = OBWs (j,k,bi,bj)
                0261      &           + obcswfac            *xx_obcsw0(j,k,bi,bj,iobcs)
                0262      &           + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs)
                0263             OBWs(j,k,bi,bj) = OBWs(j,k,bi,bj)
                0264      &           *maskW(i+ip1,j,k,bi,bj)
                0265            else if (iobcs .EQ. 3) then
                0266             OBWu(j,k,bi,bj) = OBWu (j,k,bi,bj)
                0267      &           + obcswfac            *xx_obcsw0(j,k,bi,bj,iobcs)
                0268      &           + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs)
                0269             OBWu(j,k,bi,bj) = OBWu(j,k,bi,bj)
                0270      &           *maskW(i+ip1,j,k,bi,bj)
                0271            else if (iobcs .EQ. 4) then
                0272             OBWv(j,k,bi,bj) = OBWv (j,k,bi,bj)
                0273      &           + obcswfac            *xx_obcsw0(j,k,bi,bj,iobcs)
                0274      &           + (1. _d 0 - obcswfac)*xx_obcsw1(j,k,bi,bj,iobcs)
                0275             OBWv(j,k,bi,bj) = OBWv(j,k,bi,bj)
                0276      &           *maskS(i,j,k,bi,bj)
                0277            endif
                0278           enddo
                0279          enddo
                0280         enddo
                0281        enddo
46c1d12550 Jean*0282 
7109a141b2 Patr*0283 C--   End over iobcs loop
                0284       enddo
                0285 
                0286 #endif /* ALLOW_OBCSW_CONTROL */
                0287 
46c1d12550 Jean*0288       return
7109a141b2 Patr*0289       end