#include "CTRL_CPPOPTIONS.h" #ifdef ALLOW_OBCS # include "OBCS_OPTIONS.h" #endif subroutine ctrl_getobcse( I mytime, I myiter, I mythid & ) c ================================================================== c SUBROUTINE ctrl_getobcse c ================================================================== c c o Get eastern obc of the control vector and add it c to dyn. fields c c started: heimbach@mit.edu, 29-Aug-2001 c c ================================================================== c SUBROUTINE ctrl_getobcse c ================================================================== implicit none #ifdef ALLOW_OBCSE_CONTROL c == global variables == #include "EEPARAMS.h" #include "SIZE.h" #include "PARAMS.h" #include "GRID.h" #include "OBCS.h" #include "ctrl.h" #include "ctrl_dummy.h" #include "optim.h" c == routine arguments == _RL mytime integer myiter integer mythid c == local variables == integer bi,bj integer i,j,k integer itlo,ithi integer jtlo,jthi integer jmin,jmax integer imin,imax integer ilobcse integer iobcs _RL dummy _RL obcsefac logical obcsefirst logical obcsechanged integer obcsecount0 integer obcsecount1 integer ip1 cih Integer nk,nz _RL tmpz (nr,nsx,nsy) _RL stmp character*(80) fnamein cih logical doglobalread logical ladinit character*(80) fnameobcse c == external functions == integer ilnblnk external ilnblnk c == end of interface == jtlo = mybylo(mythid) jthi = mybyhi(mythid) itlo = mybxlo(mythid) ithi = mybxhi(mythid) jmin = 1-oly jmax = sny+oly imin = 1-olx imax = snx+olx ip1 = 0 c-- Now, read the control vector. doglobalread = .false. ladinit = .false. if (optimcycle .ge. 0) then ilobcse=ilnblnk( xx_obcse_file ) write(fnameobcse(1:80),'(2a,i10.10)') & xx_obcse_file(1:ilobcse), '.', optimcycle endif c-- Get the counters, flags, and the interpolation factor. call ctrl_get_gen_rec( I xx_obcsestartdate, xx_obcseperiod, O obcsefac, obcsefirst, obcsechanged, O obcsecount0,obcsecount1, I mytime, myiter, mythid ) cih do iobcs = 1,nobcs if ( obcsefirst ) then call active_read_yz( fnameobcse, tmpfldyz, & (obcsecount0-1)*nobcs+iobcs, & doglobalread, ladinit, optimcycle, & mythid, xx_obcse_dummy ) cih if ( optimcycle .ge. 0) then cih If normal velocity if (iobcs .eq. 3) then cih Begin loop over y-points. do bj = jtlo,jthi do bi = itlo, ithi do j = jmin,jmax cih If open boundary. if ( OB_Ie(j,bi,bj) .ne. 0. ) then i = OB_Ie(j,bi,bj) cih Determine number of open vertical layers. nz = 0 do k = 1,Nr nz = nz + maskW(i+ip1,j,k,bi,bj) end do cih Compute absolute velocities from the barotropic-baroclinic modes. do k = 1,Nr if (k.le.nz) then stmp = 0. do nk = 1,nz stmp = stmp + & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj) end do tmpz(k,bi,bj) = stmp else tmpz(k,bi,bj) = 0. end if end do do k = 1,Nr tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) end do cih End open boundary. end if cih End loop over x-points end do end do end do cih End if iobcs = 3. end if cih cih If tangential velocity if (iobcs .eq. 4) then cih Begin loop over y-points. do bj = jtlo,jthi do bi = itlo, ithi do j = jmin,jmax cih If open boundary. if ( OB_Ie(j,bi,bj) .ne. 0. ) then i = OB_Ie(j,bi,bj) cih Determine number of open vertical layers. nz = 0 do k = 1,Nr nz = nz + maskS(i,j,k,bi,bj) end do cih Compute absolute velocities from the barotropic-baroclinic modes. do k = 1,Nr if (k.le.nz) then stmp = 0. do nk = 1,nz stmp = stmp + & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj) end do tmpz(k,bi,bj) = stmp else tmpz(k,bi,bj) = 0. end if end do do k = 1,Nr tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) end do cih End open boundary. end if cih End loop over x-points end do end do end do cih End if iobcs = 4. end if cih End if optimcycle > 0 . end if cih do bj = jtlo,jthi do bi = itlo,ithi do k = 1,nr do j = jmin,jmax xx_obcse1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj) enddo enddo enddo enddo cih End if obcsefirst endif cih if ( (obcsefirst) .or. (obcsechanged) ) then cih do bj = jtlo,jthi do bi = itlo,ithi do j = jmin,jmax do k = 1,nr tmpfldyz(j,k,bi,bj) = xx_obcse1(j,k,bi,bj,iobcs) enddo enddo enddo enddo cih call exf_swapffields_yz( tmpfldyz2, tmpfldyz, mythid) cih do bj = jtlo,jthi do bi = itlo,ithi do j = jmin,jmax do k = 1,nr xx_obcse0(j,k,bi,bj,iobcs) = tmpfldyz2(j,k,bi,bj) enddo enddo enddo enddo cih call active_read_yz( fnameobcse, tmpfldyz, & (obcsecount1-1)*nobcs+iobcs, & doglobalread, ladinit, optimcycle, & mythid, xx_obcse_dummy ) if ( optimcycle .ge. 0) then cih If normal velocity if (iobcs .eq. 3) then cih Begin loop over y-points. do bj = jtlo,jthi do bi = itlo, ithi do j = jmin,jmax cih If open boundary. if ( OB_Ie(j,bi,bj) .ne. 0. ) then i = OB_Ie(j,bi,bj) cih Determine number of open vertical layers. nz = 0 do k = 1,Nr nz = nz + maskW(i+ip1,j,k,bi,bj) end do cih Compute absolute velocities from the barotropic-baroclinic modes. do k = 1,Nr if (k.le.nz) then stmp = 0. do nk = 1,nz stmp = stmp + & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj) end do tmpz(k,bi,bj) = stmp else tmpz(k,bi,bj) = 0. end if end do do k = 1,Nr tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) end do cih End open boundary. end if cih End loop over x-points end do end do end do cih End if iobcs = 3. end if cih cih If tangential velocity if (iobcs .eq. 4) then cih Begin loop over y-points. do bj = jtlo,jthi do bi = itlo, ithi do j = jmin,jmax cih If open boundary. if ( OB_Ie(j,bi,bj) .ne. 0. ) then i = OB_Ie(j,bi,bj) cih Determine number of open vertical layers. nz = 0 do k = 1,Nr nz = nz + maskS(i,j,k,bi,bj) end do cih Compute absolute velocities from the barotropic-baroclinic modes. do k = 1,Nr if (k.le.nz) then stmp = 0. do nk = 1,nz stmp = stmp + & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj) end do tmpz(k,bi,bj) = stmp else tmpz(k,bi,bj) = 0. end if end do do k = 1,Nr tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj) end do cih End open boundary. end if cih End loop over x-points end do end do end do cih End if iobcs = 4. end if cih End if optimcycle > 0 . end if cih do bj = jtlo,jthi do bi = itlo,ithi do k = 1,nr do j = jmin,jmax xx_obcse1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj) enddo enddo enddo enddo cih End if obcsefirst .or. obcsechanged endif c-- Add control to model variable. do bj = jtlo,jthi do bi = itlo,ithi c-- Calculate mask for tracer cells (0 => land, 1 => water). do k = 1,nr do j = 1,sny i = OB_Ie(j,bi,bj) if (iobcs .EQ. 1) then OBEt(j,k,bi,bj) = OBEt (j,k,bi,bj) & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj) & *maskW(i+ip1,j,k,bi,bj) else if (iobcs .EQ. 2) then OBEs(j,k,bi,bj) = OBEs (j,k,bi,bj) & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj) & *maskW(i+ip1,j,k,bi,bj) else if (iobcs .EQ. 3) then OBEu(j,k,bi,bj) = OBEu (j,k,bi,bj) & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj) & *maskW(i+ip1,j,k,bi,bj) else if (iobcs .EQ. 4) then OBEv(j,k,bi,bj) = OBEv (j,k,bi,bj) & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs) & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs) OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj) & *maskS(i,j,k,bi,bj) endif enddo enddo enddo enddo C-- End over iobcs loop enddo #else /* ALLOW_OBCSE_CONTROL undefined */ c == routine arguments == _RL mytime integer myiter integer mythid c-- CPP flag ALLOW_OBCSE_CONTROL undefined. #endif /* ALLOW_OBCSE_CONTROL */ end