[3] | 1 | MODULE obc_oce |
---|
| 2 | !!============================================================================== |
---|
| 3 | !! *** MODULE obc_oce *** |
---|
| 4 | !! Open Boundary Cond. : define related variables |
---|
| 5 | !!============================================================================== |
---|
| 6 | #if defined key_obc |
---|
| 7 | !!---------------------------------------------------------------------- |
---|
| 8 | !! 'key_obc' : Open Boundary Condition |
---|
| 9 | !!---------------------------------------------------------------------- |
---|
| 10 | !! history : |
---|
| 11 | !! 8.0 01/91 (CLIPPER) Original code |
---|
| 12 | !! 8.5 06/02 (C. Talandier) modules |
---|
| 13 | !!---------------------------------------------------------------------- |
---|
| 14 | !! * Modules used |
---|
| 15 | USE par_oce ! ocean parameters |
---|
| 16 | USE obc_par ! open boundary condition parameters |
---|
| 17 | |
---|
| 18 | IMPLICIT NONE |
---|
| 19 | PUBLIC |
---|
[247] | 20 | !!---------------------------------------------------------------------- |
---|
| 21 | !! OPA 9.0 , LOCEAN-IPSL (2005) |
---|
| 22 | !! $Header$ |
---|
| 23 | !! This software is governed by the CeCILL licence see modipsl/doc/NEMO_CeCILL.txt |
---|
| 24 | !!---------------------------------------------------------------------- |
---|
[3] | 25 | |
---|
| 26 | !!---------------------------------------------------------------------- |
---|
| 27 | !! open boundary variables |
---|
| 28 | !!---------------------------------------------------------------------- |
---|
| 29 | !! |
---|
| 30 | !!General variables for open boundaries: |
---|
| 31 | !!-------------------------------------- |
---|
[32] | 32 | INTEGER :: & !: * namelist ??? * |
---|
| 33 | nbobc = 1 , & !: number of open boundaries ( 1=< nbobc =< 4 ) |
---|
| 34 | nobc_dta = 0 , & !: = 0 use the initial state as obc data |
---|
| 35 | ! ! = 1 read obc data in obcxxx.dta files |
---|
| 36 | nmoisold , & !: number of the last read month on the OBC |
---|
| 37 | nbef, naft !: index of the aftera and before fields on the OBC |
---|
[3] | 38 | |
---|
[32] | 39 | REAL(wp) :: & !!: open boundary namelist (namobc) |
---|
| 40 | rdpein = 1. , & !: damping time scale for inflow at East open boundary |
---|
| 41 | rdpwin = 1. , & !: " " at West open boundary |
---|
| 42 | rdpsin = 1. , & !: " " at South open boundary |
---|
| 43 | rdpnin = 1. , & !: " " at North open boundary |
---|
| 44 | rdpeob = 15. , & !: damping time scale for the climatology at East open boundary |
---|
| 45 | rdpwob = 15. , & !: " " at West open boundary |
---|
| 46 | rdpsob = 15. , & !: " " at South open boundary |
---|
| 47 | rdpnob = 15. , & !: " " at North open boundary |
---|
| 48 | volemp = 1. !: = 0 the total volume will have the variability of the |
---|
| 49 | ! surface Flux E-P else (volemp = 1) the volume will be constant |
---|
| 50 | ! = 1 the volume will be constant during all the integration. |
---|
[3] | 51 | |
---|
[32] | 52 | LOGICAL :: & !: |
---|
| 53 | lfbceast, lfbcwest, & !: logical flag for a fixed East and West open boundaries |
---|
| 54 | lfbcnorth, lfbcsouth !: logical flag for a fixed North and South open boundaries |
---|
| 55 | ! ! These logical flags are set to 'true' if damping time |
---|
| 56 | ! ! scale are set to 0 in the namelist, for both inflow and outflow). |
---|
[3] | 57 | |
---|
[32] | 58 | REAL(wp), DIMENSION(jpi,jpj) :: & !: |
---|
| 59 | obctmsk !: mask array identical to tmask, execpt along OBC where it is set to 0 |
---|
| 60 | ! ! it used to calculate the cumulate flux E-P in the obcvol.F90 routine |
---|
[3] | 61 | |
---|
[32] | 62 | !!---------------- |
---|
[3] | 63 | !! Rigid lid case: |
---|
| 64 | !!---------------- |
---|
[32] | 65 | INTEGER :: nbic !: number of isolated coastlines ( 0 <= nbic <= 3 ) |
---|
[3] | 66 | |
---|
[32] | 67 | INTEGER, DIMENSION(jpnic,0:4,3) :: & !: |
---|
| 68 | miic, mjic !: position of isolated coastlines points |
---|
[3] | 69 | |
---|
[32] | 70 | INTEGER, DIMENSION(0:4,3) :: & !: |
---|
| 71 | mnic !: number of points on isolated coastlines |
---|
[3] | 72 | |
---|
[32] | 73 | REAL(wp), DIMENSION(jpi,jpj) :: & !: |
---|
| 74 | gcbob !: right hand side of the barotropic elliptic equation associated |
---|
| 75 | ! ! with the OBC |
---|
[3] | 76 | |
---|
[32] | 77 | REAL(wp), DIMENSION(jpi,jpj,3) :: & !: |
---|
| 78 | gcfobc !: coef. associated with the contribution of isolated coastlines |
---|
| 79 | ! ! to the right hand side of the barotropic elliptic equation |
---|
[3] | 80 | |
---|
[32] | 81 | REAL(wp), DIMENSION(3) :: & !: |
---|
| 82 | gcbic !: time variation of the barotropic stream function along the |
---|
| 83 | ! ! isolated coastlines |
---|
[3] | 84 | |
---|
[32] | 85 | REAL(wp), DIMENSION(1) :: & !: |
---|
| 86 | bsfic0 !: barotropic stream function on isolated coastline |
---|
[3] | 87 | |
---|
[32] | 88 | REAL(wp), DIMENSION(3) :: & !: |
---|
| 89 | bsfic !: barotropic stream function on isolated coastline |
---|
[3] | 90 | |
---|
[32] | 91 | !!-------------------- |
---|
[3] | 92 | !! East open boundary: |
---|
| 93 | !!-------------------- |
---|
[32] | 94 | INTEGER :: nie0 , nie1 !: do loop index in mpp case for jpieob |
---|
| 95 | INTEGER :: nie0p1, nie1p1 !: do loop index in mpp case for jpieob+1 |
---|
| 96 | INTEGER :: nie0m1, nie1m1 !: do loop index in mpp case for jpieob-1 |
---|
| 97 | INTEGER :: nje0 , nje1 !: do loop index in mpp case for jpjed, jpjef |
---|
| 98 | INTEGER :: nje0p1, nje1m1 !: do loop index in mpp case for jpjedp1,jpjefm1 |
---|
| 99 | INTEGER :: nje1m2, nje0m1 !: do loop index in mpp case for jpjefm1-1,jpjed |
---|
[3] | 100 | |
---|
[32] | 101 | REAL(wp), DIMENSION(jpj) :: & !: |
---|
| 102 | bsfeob !: now barotropic stream fuction computed at the OBC. The corres- |
---|
| 103 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
[3] | 104 | |
---|
[32] | 105 | REAL(wp), DIMENSION(jpj,3,3) :: & !: |
---|
| 106 | bebnd !: east boundary barotropic streamfunction over 3 rows |
---|
| 107 | ! ! and 3 time step (now, before, and before before) |
---|
[3] | 108 | |
---|
[32] | 109 | REAL(wp), DIMENSION(jpjed:jpjef) :: & !: |
---|
| 110 | bfoe !: now climatology of the east boundary barotropic stream function |
---|
[3] | 111 | |
---|
[32] | 112 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
| 113 | ufoe, vfoe, & !: now climatology of the east boundary velocities |
---|
| 114 | tfoe, sfoe, & !: now climatology of the east boundary temperature and salinity |
---|
| 115 | uclie !: baroclinic componant of the zonal velocity after radiation |
---|
| 116 | ! ! in the obcdyn.F90 routine |
---|
[3] | 117 | |
---|
[32] | 118 | REAL(wp), DIMENSION(jpjglo,jpk,1) :: & !: |
---|
| 119 | uedta, tedta, sedta !: array used for interpolating monthly data on the east boundary |
---|
[3] | 120 | |
---|
[32] | 121 | !!------------------------------- |
---|
[3] | 122 | !! Arrays for radiative East OBC: |
---|
| 123 | !!------------------------------- |
---|
[32] | 124 | REAL(wp), DIMENSION(jpj,jpk,3,3) :: & !: |
---|
| 125 | uebnd, vebnd !: baroclinic u & v component of the velocity over 3 rows |
---|
[3] | 126 | ! and 3 time step (now, before, and before before) |
---|
| 127 | |
---|
[32] | 128 | REAL(wp), DIMENSION(jpj,jpk,2,2) :: & !: |
---|
| 129 | tebnd, sebnd !: East boundary temperature and salinity over 2 rows |
---|
[3] | 130 | ! and 2 time step (now and before) |
---|
| 131 | |
---|
[32] | 132 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
| 133 | u_cxebnd, v_cxebnd !: Zonal component of the phase speed ratio computed with |
---|
[3] | 134 | ! radiation of u and v velocity (respectively) at the |
---|
| 135 | ! east open boundary (u_cxebnd = cx rdt ) |
---|
| 136 | |
---|
[32] | 137 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
| 138 | uemsk, vemsk, temsk !: 2D mask for the East OB |
---|
[3] | 139 | |
---|
| 140 | ! Note that those arrays are optimized for mpp case |
---|
| 141 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
| 142 | |
---|
| 143 | !!-------------------- |
---|
[32] | 144 | !! West open boundary |
---|
| 145 | !!-------------------- |
---|
| 146 | INTEGER :: niw0 , niw1 !: do loop index in mpp case for jpiwob |
---|
| 147 | INTEGER :: niw0p1, niw1p1 !: do loop index in mpp case for jpiwob+1 |
---|
| 148 | INTEGER :: njw0 , njw1 !: do loop index in mpp case for jpjwd, jpjwf |
---|
| 149 | INTEGER :: njw0p1, njw1m1 !: do loop index in mpp case for jpjwdp1,jpjwfm1 |
---|
| 150 | INTEGER :: njw1m2, njw0m1 !: do loop index in mpp case for jpjwfm2,jpjwd |
---|
[3] | 151 | |
---|
[32] | 152 | REAL(wp), DIMENSION(jpj) :: & !: |
---|
| 153 | bsfwob !: now barotropic stream fuction computed at the OBC. The corres- |
---|
| 154 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
[3] | 155 | |
---|
[32] | 156 | REAL(wp), DIMENSION(jpj,3,3) :: & !: |
---|
| 157 | bwbnd !: West boundary barotropic streamfunction over |
---|
| 158 | ! ! 3 rows and 3 time step (now, before, and before before) |
---|
[3] | 159 | |
---|
[32] | 160 | REAL(wp), DIMENSION(jpjwd:jpjwf) :: & !: |
---|
| 161 | bfow !: now climatology of the west boundary barotropic stream function |
---|
[3] | 162 | |
---|
[32] | 163 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
| 164 | ufow, vfow, & !: now climatology of the west velocities |
---|
| 165 | tfow, sfow, & !: now climatology of the west temperature and salinity |
---|
| 166 | ucliw !: baroclinic componant of the zonal velocity after the radiation |
---|
| 167 | ! ! in the obcdyn.F90 routine |
---|
[3] | 168 | |
---|
[32] | 169 | REAL(wp), DIMENSION(jpjglo,jpk,1) :: & !: |
---|
| 170 | uwdta, twdta, swdta !: array used for interpolating monthly data on the west boundary |
---|
[3] | 171 | |
---|
| 172 | !!------------------------------- |
---|
[32] | 173 | !! Arrays for radiative West OBC |
---|
| 174 | !!------------------------------- |
---|
| 175 | REAL(wp), DIMENSION(jpj,jpk,3,3) :: & !: |
---|
| 176 | uwbnd, vwbnd !: baroclinic u & v components of the velocity over 3 rows |
---|
| 177 | ! ! and 3 time step (now, before, and before before) |
---|
[3] | 178 | |
---|
[32] | 179 | REAL(wp), DIMENSION(jpj,jpk,2,2) :: & !: |
---|
| 180 | twbnd, swbnd !: west boundary temperature and salinity over 2 rows and |
---|
| 181 | ! ! 2 time step (now and before) |
---|
[3] | 182 | |
---|
[32] | 183 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
| 184 | u_cxwbnd, v_cxwbnd !: Zonal component of the phase speed ratio computed with |
---|
| 185 | ! ! radiation of zonal and meridional velocity (respectively) |
---|
| 186 | ! ! at the west open boundary (u_cxwbnd = cx rdt ) |
---|
[3] | 187 | |
---|
[32] | 188 | REAL(wp), DIMENSION(jpj,jpk) :: & !: |
---|
| 189 | uwmsk, vwmsk, twmsk !: 2D mask for the West OB |
---|
[3] | 190 | |
---|
| 191 | ! Note that those arrays are optimized for mpp case |
---|
| 192 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
| 193 | |
---|
| 194 | !!--------------------- |
---|
[32] | 195 | !! North open boundary |
---|
| 196 | !!--------------------- |
---|
| 197 | INTEGER :: nin0 , nin1 !: do loop index in mpp case for jpind, jpinf |
---|
| 198 | INTEGER :: nin0p1, nin1m1 !: do loop index in mpp case for jpindp1, jpinfm1 |
---|
| 199 | INTEGER :: nin1m2, nin0m1 !: do loop index in mpp case for jpinfm1-1,jpind |
---|
| 200 | INTEGER :: njn0 , njn1 !: do loop index in mpp case for jpnob |
---|
| 201 | INTEGER :: njn0p1, njn1p1 !: do loop index in mpp case for jpnob+1 |
---|
| 202 | INTEGER :: njn0m1, njn1m1 !: do loop index in mpp case for jpnob-1 |
---|
[3] | 203 | |
---|
[32] | 204 | REAL(wp), DIMENSION(jpi) :: & !: |
---|
| 205 | bsfnob !: now barotropic stream fuction computed at the OBC. The corres- |
---|
| 206 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
[3] | 207 | |
---|
[32] | 208 | REAL(wp), DIMENSION(jpi,3,3) :: & !: |
---|
| 209 | bnbnd !: north boundary barotropic streamfunction over |
---|
| 210 | ! ! 3 rows and 3 time step (now, before, and before before) |
---|
[3] | 211 | |
---|
[32] | 212 | REAL(wp), DIMENSION(jpind:jpinf) :: & !: |
---|
| 213 | bfon !: now climatology of the north boundary barotropic stream function |
---|
[3] | 214 | |
---|
[32] | 215 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
| 216 | ufon, vfon, & !: now climatology of the north boundary velocities |
---|
| 217 | tfon, sfon, & !: now climatology of the north boundary temperature and salinity |
---|
| 218 | vclin !: baroclinic componant of the meridian velocity after the radiation |
---|
| 219 | ! ! in yhe obcdyn.F90 routine |
---|
[3] | 220 | |
---|
[32] | 221 | REAL(wp), DIMENSION(jpiglo,jpk,1) :: & !: |
---|
| 222 | vndta, tndta, sndta !: array used for interpolating monthly data on the north boundary |
---|
[3] | 223 | |
---|
| 224 | !!-------------------------------- |
---|
[32] | 225 | !! Arrays for radiative North OBC |
---|
| 226 | !!-------------------------------- |
---|
[3] | 227 | !! |
---|
[32] | 228 | REAL(wp), DIMENSION(jpi,jpk,3,3) :: & !: |
---|
| 229 | unbnd, vnbnd !: baroclinic u & v components of the velocity over 3 |
---|
| 230 | ! ! rows and 3 time step (now, before, and before before) |
---|
[3] | 231 | |
---|
[32] | 232 | REAL(wp), DIMENSION(jpi,jpk,2,2) :: & !: |
---|
| 233 | tnbnd, snbnd !: north boundary temperature and salinity over |
---|
| 234 | ! ! 2 rows and 2 time step (now and before) |
---|
[3] | 235 | |
---|
[32] | 236 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
| 237 | u_cynbnd, v_cynbnd !: Meridional component of the phase speed ratio compu- |
---|
| 238 | ! ! ted with radiation of zonal and meridional velocity |
---|
| 239 | ! ! (respectively) at the north OB (u_cynbnd = cx rdt ) |
---|
[3] | 240 | |
---|
[32] | 241 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
| 242 | unmsk, vnmsk, tnmsk !: 2D mask for the North OB |
---|
[3] | 243 | |
---|
| 244 | ! Note that those arrays are optimized for mpp case |
---|
| 245 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
| 246 | |
---|
| 247 | !!--------------------- |
---|
[32] | 248 | !! South open boundary |
---|
| 249 | !!--------------------- |
---|
| 250 | INTEGER :: nis0 , nis1 !: do loop index in mpp case for jpisd, jpisf |
---|
| 251 | INTEGER :: nis0p1, nis1m1 !: do loop index in mpp case for jpisdp1, jpisfm1 |
---|
| 252 | INTEGER :: nis1m2, nis0m1 !: do loop index in mpp case for jpisfm1-1,jpisd |
---|
| 253 | INTEGER :: njs0 , njs1 !: do loop index in mpp case for jpsob |
---|
| 254 | INTEGER :: njs0p1, njs1p1 !: do loop index in mpp case for jpsob+1 |
---|
[3] | 255 | |
---|
[32] | 256 | REAL(wp), DIMENSION(jpi) :: & !: |
---|
| 257 | bsfsob !: now barotropic stream fuction computed at the OBC.The corres- |
---|
| 258 | ! ! ponding bsfn will be computed by the forward time step in dynspg. |
---|
| 259 | REAL(wp), DIMENSION(jpi,3,3) :: & !: |
---|
| 260 | bsbnd !: south boundary barotropic stream function over |
---|
| 261 | ! ! 3 rows and 3 time step (now, before, and before before) |
---|
[3] | 262 | |
---|
[32] | 263 | REAL(wp), DIMENSION(jpisd:jpisf) :: & !: |
---|
| 264 | bfos !: now climatology of the south boundary barotropic stream function |
---|
[3] | 265 | |
---|
[32] | 266 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
| 267 | ufos, vfos, & !: now climatology of the south boundary velocities |
---|
| 268 | tfos, sfos, & !: now climatology of the south boundary temperature and salinity |
---|
| 269 | vclis !: baroclinic componant of the meridian velocity after the radiation |
---|
| 270 | ! ! in the obcdyn.F90 routine |
---|
[3] | 271 | |
---|
[32] | 272 | REAL(wp), DIMENSION(jpiglo,jpk,1) :: & !: |
---|
| 273 | vsdta, tsdta, ssdta !: array used for interpolating monthly data on the south boundary |
---|
[3] | 274 | |
---|
| 275 | !!-------------------------------- |
---|
[32] | 276 | !! Arrays for radiative South OBC |
---|
| 277 | !!-------------------------------- |
---|
[3] | 278 | !! computed by the forward time step in dynspg. |
---|
[32] | 279 | REAL(wp), DIMENSION(jpi,jpk,3,3) :: & !: |
---|
| 280 | usbnd, vsbnd !: baroclinic u & v components of the velocity over 3 |
---|
| 281 | ! ! rows and 3 time step (now, before, and before before) |
---|
[3] | 282 | |
---|
[32] | 283 | REAL(wp), DIMENSION(jpi,jpk,2,2) :: & !: |
---|
| 284 | tsbnd, ssbnd !: south boundary temperature and salinity over |
---|
| 285 | ! ! 2 rows and 2 time step (now and before) |
---|
[3] | 286 | |
---|
[32] | 287 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
| 288 | u_cysbnd, v_cysbnd !: Meridional component of the phase speed ratio compu- |
---|
| 289 | ! ! ted with radiation of zonal and meridional velocity |
---|
| 290 | ! ! (repsectively) at the south OB (u_cynbnd = cx rdt ) |
---|
[3] | 291 | |
---|
[32] | 292 | REAL(wp), DIMENSION(jpi,jpk) :: & !: |
---|
| 293 | usmsk, vsmsk, tsmsk !: 2D mask for the South OB |
---|
[3] | 294 | |
---|
| 295 | ! Note that those arrays are optimized for mpp case |
---|
| 296 | ! (hence the dimension jpj is the size of one processor subdomain) |
---|
| 297 | |
---|
| 298 | #else |
---|
| 299 | !!---------------------------------------------------------------------- |
---|
| 300 | !! Default option : Empty module |
---|
| 301 | !!---------------------------------------------------------------------- |
---|
| 302 | #endif |
---|
[32] | 303 | |
---|
[3] | 304 | !!====================================================================== |
---|
| 305 | END MODULE obc_oce |
---|