File indexing completed on 2024-06-06 05:10:34 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_getobcss(
9f5240b52a Jean*0007 I myTime, myIter, myThid )
7109a141b2 Patr*0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020
0021
0022
0023 implicit none
0024
0025
46c1d12550 Jean*0026 #ifdef ALLOW_OBCSS_CONTROL
7109a141b2 Patr*0027 #include "EEPARAMS.h"
0028 #include "SIZE.h"
0029 #include "PARAMS.h"
0030 #include "GRID.h"
af61e5eb16 Mart*0031 #include "OBCS_PARAMS.h"
46c1d12550 Jean*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_OBCSS_CONTROL */
7109a141b2 Patr*0040
0041
9f5240b52a Jean*0042 _RL myTime
0043 integer myIter
0044 integer myThid
7109a141b2 Patr*0045
46c1d12550 Jean*0046 #ifdef ALLOW_OBCSS_CONTROL
9f5240b52a Jean*0047
0048 integer ilnblnk
0049 external ilnblnk
7109a141b2 Patr*0050
9f5240b52a Jean*0051
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 ilobcss
0058 integer iobcs
0059 _RL obcssfac
0060 logical obcssfirst
0061 logical obcsschanged
0062 integer obcsscount0
0063 integer obcsscount1
0064 integer jp1
9f5240b52a Jean*0065
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) fnameobcss
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
0078
9f5240b52a Jean*0079 jtlo = myByLo(myThid)
0080 jthi = myByHi(myThid)
0081 itlo = myBxLo(myThid)
0082 ithi = myBxHi(myThid)
0083 imin = 1-OLx
0084 imax = sNx+OLx
7109a141b2 Patr*0085 jp1 = 1
0086
0087
0088 doglobalread = .false.
0089 ladinit = .false.
0090
f9d7cbfb72 Ou W*0091
0092 ilDir = ilnblnk(ctrlDir)
0093
7109a141b2 Patr*0094 if (optimcycle .ge. 0) then
859d62c7fb Mart*0095 ilobcss=ilnblnk( xx_obcss_file )
de57a2ec4b Mart*0096 write(fnameobcss,'(2a,i10.10)')
f9d7cbfb72 Ou W*0097 & ctrlDir(1:ilDir)//xx_obcss_file(1:ilobcss), '.', optimcycle
7109a141b2 Patr*0098 endif
0099
0100
0101 call ctrl_get_gen_rec(
0102 I xx_obcssstartdate, xx_obcssperiod,
0103 O obcssfac, obcssfirst, obcsschanged,
0104 O obcsscount0,obcsscount1,
9f5240b52a Jean*0105 I myTime, myIter, myThid )
7109a141b2 Patr*0106
0107 do iobcs = 1,nobcs
9f5240b52a Jean*0108
859d62c7fb Mart*0109 if ( obcssfirst ) then
720a211d38 Ou W*0110 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0111 call active_read_xz( fnameobcss, tmpfldxz,
0112 & (obcsscount0-1)*nobcs+iobcs,
0113 & doglobalread, ladinit, optimcycle,
4d72283393 Mart*0114 & myThid, xx_obcss_dummy )
720a211d38 Ou W*0115 #else
0116 CALL READ_REC_XZ_RL( fnameobcss, ctrlprec, Nr, tmpfldxz,
0117 & (obcsscount0-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 i = imin,imax
0125 j = OB_Js(i,bi,bj)
9fd29e65a3 Jean*0126 IF ( j.EQ.OB_indexNone ) j = 1
b6199164d9 Matt*0127
0128 nz = 0
0129 do k = 1,Nr
0130 if (iobcs .eq. 3) then
0131 nz = nz + maskS(i,j+jp1,k,bi,bj)
0132 else
0133 nz = nz + maskW(i,j,k,bi,bj)
0134 endif
0135 end do
0136
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)*tmpfldxz(i,nk,bi,bj)
0143 end do
0144 tmpz(k,bi,bj) = stmp
0145 else
0146 tmpz(k,bi,bj) = 0.
0147 end if
0148 end do
0149 do k = 1,Nr
0150 if (iobcs .eq. 3) then
0151 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
0152 & *recip_hFacS(i,j+jp1,k,bi,bj)
0153 else
0154 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
0155 & *recip_hFacW(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 i = imin,imax
0163 xx_obcss1(i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
0164
7c7f2df94e Mart*0165 enddo
7109a141b2 Patr*0166 enddo
859d62c7fb Mart*0167 enddo
0168 enddo
0169 endif
46c1d12550 Jean*0170
859d62c7fb Mart*0171 if ( (obcssfirst) .or. (obcsschanged)) 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 i = imin,imax
0177 xx_obcss0(i,k,bi,bj,iobcs) = xx_obcss1(i,k,bi,bj,iobcs)
0178 tmpfldxz (i,k,bi,bj) = 0. _d 0
0179 enddo
7109a141b2 Patr*0180 enddo
859d62c7fb Mart*0181 enddo
0182 enddo
0183
720a211d38 Ou W*0184 #ifdef ALLOW_AUTODIFF
859d62c7fb Mart*0185 call active_read_xz( fnameobcss, tmpfldxz,
0186 & (obcsscount1-1)*nobcs+iobcs,
0187 & doglobalread, ladinit, optimcycle,
4d72283393 Mart*0188 & myThid, xx_obcss_dummy )
720a211d38 Ou W*0189 #else
0190 CALL READ_REC_XZ_RL( fnameobcss, ctrlprec, Nr, tmpfldxz,
0191 & (obcsscount1-1)*nobcs+iobcs, 1, myThid )
0192 #endif
7109a141b2 Patr*0193
0194 do bj = jtlo,jthi
7aa90384e1 Mart*0195 do bi = itlo,ithi
b6199164d9 Matt*0196 #ifdef ALLOW_OBCS_CONTROL_MODES
0197 if (iobcs .gt. 2) then
0198 do i = imin,imax
0199 j = OB_Js(i,bi,bj)
9fd29e65a3 Jean*0200 IF ( j.EQ.OB_indexNone ) j = 1
b6199164d9 Matt*0201
0202 nz = 0
0203 do k = 1,Nr
0204 if (iobcs .eq. 3) then
0205 nz = nz + maskS(i,j+jp1,k,bi,bj)
0206 else
0207 nz = nz + maskW(i,j,k,bi,bj)
0208 endif
0209 end do
0210
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)*tmpfldxz(i,nk,bi,bj)
0217 end do
0218 tmpz(k,bi,bj) = stmp
0219 else
0220 tmpz(k,bi,bj) = 0.
0221 end if
0222 end do
0223 do k = 1,Nr
0224 if (iobcs .eq. 3) then
0225 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
0226 & *recip_hFacS(i,j+jp1,k,bi,bj)
0227 else
0228 tmpfldxz(i,k,bi,bj) = tmpz(k,bi,bj)
0229 & *recip_hFacW(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 i = imin,imax
0237 xx_obcss1 (i,k,bi,bj,iobcs) = tmpfldxz (i,k,bi,bj)
0238
7109a141b2 Patr*0239 enddo
7aa90384e1 Mart*0240 enddo
0241 enddo
7109a141b2 Patr*0242 enddo
859d62c7fb Mart*0243 endif
0244
0245
0246 do bj = jtlo,jthi
0247 do bi = itlo,ithi
0248
9f5240b52a Jean*0249 do k = 1,Nr
0250 do i = 1,sNx
0251 j = OB_Js(i,bi,bj)
9fd29e65a3 Jean*0252 IF ( j.EQ.OB_indexNone ) j = 1
859d62c7fb Mart*0253 if (iobcs .EQ. 1) then
0254 OBSt(i,k,bi,bj) = OBSt (i,k,bi,bj)
0255 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
0256 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
0257 OBSt(i,k,bi,bj) = OBSt(i,k,bi,bj)
0258 & *maskS(i,j+jp1,k,bi,bj)
0259 else if (iobcs .EQ. 2) then
0260 OBSs(i,k,bi,bj) = OBSs (i,k,bi,bj)
0261 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
0262 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
0263 OBSs(i,k,bi,bj) = OBSs(i,k,bi,bj)
0264 & *maskS(i,j+jp1,k,bi,bj)
0265 else if (iobcs .EQ. 4) then
0266 OBSu(i,k,bi,bj) = OBSu (i,k,bi,bj)
0267 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
0268 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
0269 OBSu(i,k,bi,bj) = OBSu(i,k,bi,bj)
0270 & *maskW(i,j,k,bi,bj)
0271 else if (iobcs .EQ. 3) then
0272 OBSv(i,k,bi,bj) = OBSv (i,k,bi,bj)
0273 & + obcssfac *xx_obcss0(i,k,bi,bj,iobcs)
0274 & + (1. _d 0 - obcssfac)*xx_obcss1(i,k,bi,bj,iobcs)
0275 OBSv(i,k,bi,bj) = OBSv(i,k,bi,bj)
0276 & *maskS(i,j+jp1,k,bi,bj)
0277 endif
0278 enddo
0279 enddo
0280 enddo
0281 enddo
7109a141b2 Patr*0282
0283
0284 enddo
0285
0286 #endif /* ALLOW_OBCSS_CONTROL */
0287
46c1d12550 Jean*0288 return
7109a141b2 Patr*0289 end