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