File indexing completed on 2024-06-06 05:10:33 UTC
view on githubraw file Latest commit af61e5eb on 2024-06-06 03:30:35 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_getobcse(
9f5240b52a Jean*0007 I myTime, myIter, myThid )
7109a141b2 Patr*0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022 implicit none
0023
0024
46c1d12550 Jean*0025 #ifdef ALLOW_OBCSE_CONTROL
7109a141b2 Patr*0026 #include "EEPARAMS.h"
0027 #include "SIZE.h"
0028 #include "PARAMS.h"
0029 #include "GRID.h"
af61e5eb16 Mart*0030 #include "OBCS_PARAMS.h"
46c1d12550 Jean*0031 #include "OBCS_GRID.h"
0032 #include "OBCS_FIELDS.h"
c04085ad02 Patr*0033 #include "CTRL_SIZE.h"
4d72283393 Mart*0034 #include "CTRL.h"
edcd27be69 Mart*0035 #include "CTRL_DUMMY.h"
65754df434 Mart*0036 #include "OPTIMCYCLE.h"
e612621177 Gael*0037 #include "CTRL_OBCS.h"
46c1d12550 Jean*0038 #endif /* ALLOW_OBCSE_CONTROL */
7109a141b2 Patr*0039
0040
9f5240b52a Jean*0041 _RL myTime
0042 integer myIter
0043 integer myThid
7109a141b2 Patr*0044
46c1d12550 Jean*0045 #ifdef ALLOW_OBCSE_CONTROL
9f5240b52a Jean*0046
0047 integer ilnblnk
0048 external ilnblnk
7109a141b2 Patr*0049
9f5240b52a Jean*0050
7109a141b2 Patr*0051 integer bi,bj
0052 integer i,j,k
0053 integer itlo,ithi
0054 integer jtlo,jthi
0055 integer jmin,jmax
0056 integer ilobcse
0057 integer iobcs
0058 _RL obcsefac
0059 logical obcsefirst
0060 logical obcsechanged
0061 integer obcsecount0
0062 integer obcsecount1
0063 integer ip1
9f5240b52a Jean*0064
0065 _RL tmpfldyz (1-OLy:sNy+OLy,Nr,nSx,nSy)
7109a141b2 Patr*0066 logical doglobalread
0067 logical ladinit
de57a2ec4b Mart*0068 character*(MAX_LEN_FNAM) fnameobcse
b6199164d9 Matt*0069 #ifdef ALLOW_OBCS_CONTROL_MODES
0070 integer nk,nz
9f5240b52a Jean*0071 _RL tmpz (Nr,nSx,nSy)
b6199164d9 Matt*0072 _RL stmp
0073 #endif
f9d7cbfb72 Ou W*0074 integer ilDir
7109a141b2 Patr*0075
0076
0077
9f5240b52a Jean*0078 jtlo = myByLo(myThid)
0079 jthi = myByHi(myThid)
0080 itlo = myBxLo(myThid)
0081 ithi = myBxHi(myThid)
0082 jmin = 1-OLy
0083 jmax = sNy+OLy
7109a141b2 Patr*0084 ip1 = 0
0085
0086
0087 doglobalread = .false.
0088 ladinit = .false.
0089
f9d7cbfb72 Ou W*0090
0091 ilDir = ilnblnk(ctrlDir)
0092
7109a141b2 Patr*0093 if (optimcycle .ge. 0) then
859d62c7fb Mart*0094 ilobcse=ilnblnk( xx_obcse_file )
de57a2ec4b Mart*0095 write(fnameobcse,'(2a,i10.10)')
f9d7cbfb72 Ou W*0096 & ctrlDir(1:ilDir)//xx_obcse_file(1:ilobcse), '.', optimcycle
7109a141b2 Patr*0097 endif
0098
0099
0100 call ctrl_get_gen_rec(
0101 I xx_obcsestartdate, xx_obcseperiod,
0102 O obcsefac, obcsefirst, obcsechanged,
0103 O obcsecount0,obcsecount1,
9f5240b52a Jean*0104 I myTime, myIter, myThid )
7109a141b2 Patr*0105
0106 do iobcs = 1,nobcs
0107
859d62c7fb Mart*0108 if ( obcsefirst ) then
720a211d38 Ou W*0109 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0110 call active_read_yz( fnameobcse, tmpfldyz,
0111 & (obcsecount0-1)*nobcs+iobcs,
0112 & doglobalread, ladinit, optimcycle,
4d72283393 Mart*0113 & myThid, xx_obcse_dummy )
720a211d38 Ou W*0114 #else
0115 CALL READ_REC_YZ_RL( fnameobcse, ctrlprec, Nr, tmpfldyz,
0116 & (obcsecount0-1)*nobcs+iobcs, 1, myThid )
0117 #endif
859d62c7fb Mart*0118
0119 do bj = jtlo,jthi
0120 do bi = itlo,ithi
b6199164d9 Matt*0121 #ifdef ALLOW_OBCS_CONTROL_MODES
0122 if (iobcs .gt. 2) then
0123 do j = jmin,jmax
0124 i = OB_Ie(j,bi,bj)
9fd29e65a3 Jean*0125 IF ( i.EQ.OB_indexNone ) i = 1
b6199164d9 Matt*0126
0127 nz = 0
0128 do k = 1,Nr
0129 if (iobcs .eq. 3) then
0130 nz = nz + maskW(i+ip1,j,k,bi,bj)
0131 else
0132 nz = nz + maskS(i,j,k,bi,bj)
0133 endif
0134 end do
0135
0136 do k = 1,Nr
0137 if (k.le.nz) then
0138 stmp = 0.
0139 do nk = 1,nz
0140 stmp = stmp +
0141 & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj)
0142 end do
0143 tmpz(k,bi,bj) = stmp
0144 else
0145 tmpz(k,bi,bj) = 0.
0146 end if
0147 enddo
0148 do k = 1,Nr
0149 if (iobcs .eq. 3) then
0150 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
0151 & *recip_hFacW(i+ip1,j,k,bi,bj)
0152 else
0153 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
0154 & *recip_hFacS(i,j,k,bi,bj)
0155 endif
0156 end do
0157 enddo
0158 endif
0159 #endif
9f5240b52a Jean*0160 do k = 1,Nr
859d62c7fb Mart*0161 do j = jmin,jmax
0162 xx_obcse1(j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
0163
7c7f2df94e Mart*0164 enddo
7109a141b2 Patr*0165 enddo
859d62c7fb Mart*0166 enddo
0167 enddo
0168 endif
46c1d12550 Jean*0169
859d62c7fb Mart*0170 if ( (obcsefirst) .or. (obcsechanged)) then
46c1d12550 Jean*0171
859d62c7fb Mart*0172 do bj = jtlo,jthi
0173 do bi = itlo,ithi
0174 do j = jmin,jmax
9f5240b52a Jean*0175 do k = 1,Nr
859d62c7fb Mart*0176 xx_obcse0(j,k,bi,bj,iobcs) = xx_obcse1(j,k,bi,bj,iobcs)
0177 tmpfldyz (j,k,bi,bj) = 0. _d 0
0178 enddo
7109a141b2 Patr*0179 enddo
859d62c7fb Mart*0180 enddo
0181 enddo
46c1d12550 Jean*0182
720a211d38 Ou W*0183 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0184 call active_read_yz( fnameobcse, tmpfldyz,
0185 & (obcsecount1-1)*nobcs+iobcs,
0186 & doglobalread, ladinit, optimcycle,
4d72283393 Mart*0187 & myThid, xx_obcse_dummy )
720a211d38 Ou W*0188 #else
0189 CALL READ_REC_YZ_RL( fnameobcse, ctrlprec, Nr, tmpfldyz,
0190 & (obcsecount1-1)*nobcs+iobcs, 1, myThid )
0191 #endif
7109a141b2 Patr*0192
0193 do bj = jtlo,jthi
7aa90384e1 Mart*0194 do bi = itlo,ithi
b6199164d9 Matt*0195 #ifdef ALLOW_OBCS_CONTROL_MODES
0196 if (iobcs .gt. 2) then
0197 do j = jmin,jmax
0198 i = OB_Ie(j,bi,bj)
9fd29e65a3 Jean*0199 IF ( i.EQ.OB_indexNone ) i = 1
b6199164d9 Matt*0200
0201 nz = 0
0202 do k = 1,Nr
0203 if (iobcs .eq. 3) then
0204 nz = nz + maskW(i+ip1,j,k,bi,bj)
0205 else
0206 nz = nz + maskS(i,j,k,bi,bj)
0207 endif
0208 end do
0209
0210 do k = 1,Nr
0211 if (k.le.nz) then
0212 stmp = 0.
0213 do nk = 1,nz
0214 stmp = stmp +
0215 & modesv(k,nk,nz)*tmpfldyz(j,nk,bi,bj)
0216 end do
0217 tmpz(k,bi,bj) = stmp
0218 else
0219 tmpz(k,bi,bj) = 0.
0220 endif
0221 enddo
0222 do k = 1,Nr
0223 if (iobcs .eq. 3) then
0224 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
0225 & *recip_hFacW(i+ip1,j,k,bi,bj)
0226 else
0227 tmpfldyz(j,k,bi,bj) = tmpz(k,bi,bj)
0228 & *recip_hFacS(i,j,k,bi,bj)
0229 endif
0230 end do
0231 enddo
0232 endif
0233 #endif
9f5240b52a Jean*0234 do k = 1,Nr
859d62c7fb Mart*0235 do j = jmin,jmax
0236 xx_obcse1 (j,k,bi,bj,iobcs) = tmpfldyz (j,k,bi,bj)
0237
7109a141b2 Patr*0238 enddo
7aa90384e1 Mart*0239 enddo
0240 enddo
7109a141b2 Patr*0241 enddo
859d62c7fb Mart*0242 endif
46c1d12550 Jean*0243
859d62c7fb Mart*0244
0245 do bj = jtlo,jthi
0246 do bi = itlo,ithi
0247
9f5240b52a Jean*0248 do k = 1,Nr
0249 do j = 1,sNy
859d62c7fb Mart*0250 i = OB_Ie(j,bi,bj)
9fd29e65a3 Jean*0251 IF ( i.EQ.OB_indexNone ) i = 1
859d62c7fb Mart*0252 if (iobcs .EQ. 1) then
0253 OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
0254 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
0255 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
0256 OBEt(j,k,bi,bj) = OBEt(j,k,bi,bj)
0257 & *maskW(i+ip1,j,k,bi,bj)
0258 else if (iobcs .EQ. 2) then
0259 OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
0260 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
0261 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
0262 OBEs(j,k,bi,bj) = OBEs(j,k,bi,bj)
0263 & *maskW(i+ip1,j,k,bi,bj)
0264 else if (iobcs .EQ. 3) then
0265 OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
0266 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
0267 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
0268 OBEu(j,k,bi,bj) = OBEu(j,k,bi,bj)
0269 & *maskW(i+ip1,j,k,bi,bj)
0270 else if (iobcs .EQ. 4) then
0271 OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
0272 & + obcsefac *xx_obcse0(j,k,bi,bj,iobcs)
0273 & + (1. _d 0 - obcsefac)*xx_obcse1(j,k,bi,bj,iobcs)
0274 OBEv(j,k,bi,bj) = OBEv(j,k,bi,bj)
0275 & *maskS(i,j,k,bi,bj)
0276 endif
0277 enddo
0278 enddo
0279 enddo
0280 enddo
46c1d12550 Jean*0281
7109a141b2 Patr*0282
0283 enddo
0284
0285 #endif /* ALLOW_OBCSE_CONTROL */
0286
46c1d12550 Jean*0287 return
7109a141b2 Patr*0288 end