#include "COST_CPPOPTIONS.h" #ifdef ALLOW_OBCS # include "OBCS_OPTIONS.h" #endif subroutine cost_obcsn( I myiter, I mytime, I startrec, I endrec, I mythid & ) c ================================================================== c SUBROUTINE cost_obcsn c ================================================================== c c o cost function contribution obc c c o G. Gebbie, gebbie@mit.edu, 18-Mar-2003 c ================================================================== c SUBROUTINE cost_obcsn c ================================================================== implicit none c == global variables == #include "EEPARAMS.h" #include "SIZE.h" #include "PARAMS.h" #include "GRID.h" #include "DYNVARS.h" #ifdef ALLOW_OBCS # include "OBCS.h" #endif #include "cal.h" #include "ecco_cost.h" #include "ctrl.h" #include "ctrl_dummy.h" #include "optim.h" c == routine arguments == integer myiter _RL mytime integer mythid cgg( integer startrec integer endrec cgg) c == local variables == integer bi,bj integer i,j,k integer itlo,ithi integer jtlo,jthi integer jmin,jmax integer imin,imax integer irec integer il integer iobcs integer jp1 integer nrec integer ilfld integer igg integer ihh _RL fctile _RL fcthread _RL dummy _RL gg _RL tmpx _RL tmpfield (1-olx:snx+olx,nr,nsx,nsy) _RL maskxz (1-olx:snx+olx,nr,nsx,nsy) _RL area _RL volflux character*(80) fnamefld logical doglobalread logical ladinit #ifdef ECCO_VERBOSE character*(MAX_LEN_MBUF) msgbuf #endif 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 jmax = sny imin = 1 imax = snx c-- Read tiled data. doglobalread = .false. ladinit = .false. c Number of records to be used. nrec = endrec-startrec+1 #ifdef ALLOW_OBCSN_COST_CONTRIBUTION jp1 = 0 fcthread = 0. _d 0 #ifdef ECCO_VERBOSE _BEGIN_MASTER( mythid ) write(msgbuf,'(a)') ' ' call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) write(msgbuf,'(a)') ' ' call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) write(msgbuf,'(a,i9.8)') & ' cost_obcsn: number of records to process: ',nrec call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) write(msgbuf,'(a)') ' ' call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) _END_MASTER( mythid ) #endif if (optimcycle .ge. 0) then ilfld=ilnblnk( xx_obcsn_file ) write(fnamefld(1:80),'(2a,i10.10)') & xx_obcsn_file(1:ilfld), '.', optimcycle endif c-- Loop over records. do irec = 1,nrec area = 0. _d 0 volflux = 0. _d 0 call active_read_xz( fnamefld, tmpfield, irec, doglobalread, & ladinit, optimcycle, mythid & , xx_obcsn_dummy ) cgg Need to solve for iobcs would have been. gg = (irec-1)/nobcs igg = int(gg) iobcs = irec - igg*nobcs call active_read_xz( 'maskobcsn', maskxz, & iobcs, & doglobalread, ladinit, 0, & mythid, dummy ) #ifdef BALANCE_CONTROL_VOLFLUX_GLOBAL cih -- Balance net transport from the northern boundary. c Compute total net transport. if (iobcs .eq. 3) then ihh = igg+1 call ctrl_volflux( ihh, area, volflux, mythid) _GLOBAL_SUM_RL( volflux, mythid ) _GLOBAL_SUM_RL( area,mythid ) c print*,'volflux,area',volflux,area endif c Correct barofield if normal velocity at the northern boundary. do bj = jtlo,jthi do bi = itlo,ithi do i = imin,imax if (iobcs .eq. 3) then tmpfield(i,1,bi,bj) = (tmpfield(i,1,bi,bj) + & volflux/area)*maskxz(i,1,bi,bj) c print*,'volflux2,area2',volflux,area c print*,'barofield',tmpfield(i,1,bi,bj) endif enddo enddo enddo #endif c-- Loop over this thread's tiles. do bj = jtlo,jthi do bi = itlo,ithi c-- Determine the weights to be used. fctile = 0. _d 0 do k = 1, Nr do i = imin,imax j = OB_Jn(I,bi,bj) cgg if (maskS(i,j+jp1,k,bi,bj) .ne. 0.) then tmpx = tmpfield(i,k,bi,bj) CMM fctile = fctile + wobcsn2(i,k,bi,bj,iobcs) fctile = fctile + wobcsn(k,iobcs) & *tmpx*tmpx*maskxz(i,k,bi,bj) cgg endif CMM if (wobcsn2(i,k,bi,bj,iobcs)*maskxz(i,k,bi,bj).ne.0.) if (wobcsn(k,iobcs)*maskxz(i,k,bi,bj).ne.0.) & num_obcsn(bi,bj) = num_obcsn(bi,bj) + 1. _d 0 cgg print*,'S fctile',fctile enddo enddo objf_obcsn(bi,bj) = objf_obcsn(bi,bj) + fctile fcthread = fcthread + fctile enddo enddo #ifdef ECCO_VERBOSE c-- Print cost function for all tiles. _GLOBAL_SUM_RL( fcthread , myThid ) write(msgbuf,'(a)') ' ' call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) write(msgbuf,'(a,i8.8)') & ' cost_obcsn: irec = ',irec call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) write(msgbuf,'(a,a,d22.15)') & ' global cost function value', & ' (obcsn) = ',fcthread call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) write(msgbuf,'(a)') ' ' call print_message( msgbuf, standardmessageunit, & SQUEEZE_RIGHT , mythid) #endif enddo c-- End of loop over records. #endif return end