[1682] | 1 | !> @file surface_coupler.f90 |
---|
[2000] | 2 | !------------------------------------------------------------------------------! |
---|
[1036] | 3 | ! This file is part of PALM. |
---|
| 4 | ! |
---|
[2000] | 5 | ! PALM is free software: you can redistribute it and/or modify it under the |
---|
| 6 | ! terms of the GNU General Public License as published by the Free Software |
---|
| 7 | ! Foundation, either version 3 of the License, or (at your option) any later |
---|
| 8 | ! version. |
---|
[1036] | 9 | ! |
---|
| 10 | ! PALM is distributed in the hope that it will be useful, but WITHOUT ANY |
---|
| 11 | ! WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR |
---|
| 12 | ! A PARTICULAR PURPOSE. See the GNU General Public License for more details. |
---|
| 13 | ! |
---|
| 14 | ! You should have received a copy of the GNU General Public License along with |
---|
| 15 | ! PALM. If not, see <http://www.gnu.org/licenses/>. |
---|
| 16 | ! |
---|
[1818] | 17 | ! Copyright 1997-2016 Leibniz Universitaet Hannover |
---|
[2000] | 18 | !------------------------------------------------------------------------------! |
---|
[1036] | 19 | ! |
---|
[258] | 20 | ! Current revisions: |
---|
[1092] | 21 | ! ------------------ |
---|
[1321] | 22 | ! |
---|
[2032] | 23 | ! |
---|
[1321] | 24 | ! Former revisions: |
---|
| 25 | ! ----------------- |
---|
| 26 | ! $Id: surface_coupler.f90 2032 2016-10-21 15:13:51Z maronga $ |
---|
| 27 | ! |
---|
[2032] | 28 | ! 2031 2016-10-21 15:11:58Z knoop |
---|
| 29 | ! renamed variable rho to rho_ocean |
---|
| 30 | ! |
---|
[2001] | 31 | ! 2000 2016-08-20 18:09:15Z knoop |
---|
| 32 | ! Forced header and separation lines into 80 columns |
---|
| 33 | ! |
---|
[1683] | 34 | ! 1682 2015-10-07 23:56:08Z knoop |
---|
| 35 | ! Code annotations made doxygen readable |
---|
| 36 | ! |
---|
[1428] | 37 | ! 1427 2014-07-07 14:04:59Z maronga |
---|
| 38 | ! Bugfix: value of l_v corrected. |
---|
| 39 | ! |
---|
[1419] | 40 | ! 1418 2014-06-06 13:05:08Z fricke |
---|
| 41 | ! Bugfix: For caluclation of the salinity flux at the sea surface, |
---|
| 42 | ! the given value for salinity must be in percent and not in psu |
---|
| 43 | ! |
---|
[1354] | 44 | ! 1353 2014-04-08 15:21:23Z heinze |
---|
| 45 | ! REAL constants provided with KIND-attribute |
---|
| 46 | ! |
---|
[1325] | 47 | ! 1324 2014-03-21 09:13:16Z suehring |
---|
| 48 | ! Bugfix: ONLY statement for module pegrid removed |
---|
| 49 | ! |
---|
[1323] | 50 | ! 1322 2014-03-20 16:38:49Z raasch |
---|
| 51 | ! REAL constants defined as wp-kind |
---|
| 52 | ! |
---|
[1321] | 53 | ! 1320 2014-03-20 08:40:49Z raasch |
---|
[1320] | 54 | ! ONLY-attribute added to USE-statements, |
---|
| 55 | ! kind-parameters added to all INTEGER and REAL declaration statements, |
---|
| 56 | ! kinds are defined in new module kinds, |
---|
| 57 | ! old module precision_kind is removed, |
---|
| 58 | ! revision history before 2012 removed, |
---|
| 59 | ! comment fields (!:) to be used for variable explanations added to |
---|
| 60 | ! all variable declaration statements |
---|
[102] | 61 | ! |
---|
[1319] | 62 | ! 1318 2014-03-17 13:35:16Z raasch |
---|
| 63 | ! module interfaces removed |
---|
| 64 | ! |
---|
[1093] | 65 | ! 1092 2013-02-02 11:24:22Z raasch |
---|
| 66 | ! unused variables removed |
---|
| 67 | ! |
---|
[1037] | 68 | ! 1036 2012-10-22 13:43:42Z raasch |
---|
| 69 | ! code put under GPL (PALM 3.9) |
---|
| 70 | ! |
---|
[881] | 71 | ! 880 2012-04-13 06:28:59Z raasch |
---|
| 72 | ! Bugfix: preprocessor statements for parallel execution added |
---|
| 73 | ! |
---|
[110] | 74 | ! 109 2007-08-28 15:26:47Z letzel |
---|
[102] | 75 | ! Initial revision |
---|
| 76 | ! |
---|
| 77 | ! Description: |
---|
| 78 | ! ------------ |
---|
[1682] | 79 | !> Data exchange at the interface between coupled models |
---|
[102] | 80 | !------------------------------------------------------------------------------! |
---|
[1682] | 81 | SUBROUTINE surface_coupler |
---|
| 82 | |
---|
[102] | 83 | |
---|
[1320] | 84 | USE arrays_3d, & |
---|
[2031] | 85 | ONLY: pt, shf, qsws, qswst_remote, rho_ocean, sa, saswst, total_2d_a, & |
---|
[1320] | 86 | total_2d_o, tswst, u, usws, uswst, v, vsws, vswst |
---|
| 87 | |
---|
[1427] | 88 | USE cloud_parameters, & |
---|
| 89 | ONLY: cp, l_v |
---|
| 90 | |
---|
[1320] | 91 | USE control_parameters, & |
---|
| 92 | ONLY: coupling_mode, coupling_mode_remote, coupling_topology, & |
---|
| 93 | humidity, humidity_remote, message_string, terminate_coupled, & |
---|
| 94 | terminate_coupled_remote, time_since_reference_point |
---|
| 95 | |
---|
| 96 | USE cpulog, & |
---|
| 97 | ONLY: cpu_log, log_point |
---|
| 98 | |
---|
| 99 | USE indices, & |
---|
| 100 | ONLY: nbgp, nx, nxl, nxlg, nxr, nxrg, nx_a, nx_o, ny, nyn, nyng, nys, & |
---|
| 101 | nysg, ny_a, ny_o, nzt |
---|
| 102 | |
---|
| 103 | USE kinds |
---|
| 104 | |
---|
[102] | 105 | USE pegrid |
---|
| 106 | |
---|
| 107 | IMPLICIT NONE |
---|
| 108 | |
---|
[1682] | 109 | REAL(wp) :: time_since_reference_point_rem !< |
---|
| 110 | REAL(wp) :: total_2d(-nbgp:ny+nbgp,-nbgp:nx+nbgp) !< |
---|
[102] | 111 | |
---|
[1682] | 112 | REAL(wp) :: cpw = 4218.0_wp !< heat capacity of water at constant pressure |
---|
[1427] | 113 | |
---|
[206] | 114 | #if defined( __parallel ) |
---|
[102] | 115 | |
---|
[667] | 116 | CALL cpu_log( log_point(39), 'surface_coupler', 'start' ) |
---|
[102] | 117 | |
---|
[667] | 118 | |
---|
| 119 | |
---|
[102] | 120 | ! |
---|
[108] | 121 | !-- In case of model termination initiated by the remote model |
---|
| 122 | !-- (terminate_coupled_remote > 0), initiate termination of the local model. |
---|
| 123 | !-- The rest of the coupler must then be skipped because it would cause an MPI |
---|
| 124 | !-- intercomminucation hang. |
---|
| 125 | !-- If necessary, the coupler will be called at the beginning of the next |
---|
| 126 | !-- restart run. |
---|
[667] | 127 | |
---|
| 128 | IF ( coupling_topology == 0 ) THEN |
---|
[709] | 129 | CALL MPI_SENDRECV( terminate_coupled, 1, MPI_INTEGER, target_id, & |
---|
| 130 | 0, & |
---|
| 131 | terminate_coupled_remote, 1, MPI_INTEGER, target_id, & |
---|
[667] | 132 | 0, comm_inter, status, ierr ) |
---|
| 133 | ELSE |
---|
| 134 | IF ( myid == 0) THEN |
---|
| 135 | CALL MPI_SENDRECV( terminate_coupled, 1, MPI_INTEGER, & |
---|
| 136 | target_id, 0, & |
---|
| 137 | terminate_coupled_remote, 1, MPI_INTEGER, & |
---|
| 138 | target_id, 0, & |
---|
| 139 | comm_inter, status, ierr ) |
---|
| 140 | ENDIF |
---|
[709] | 141 | CALL MPI_BCAST( terminate_coupled_remote, 1, MPI_INTEGER, 0, comm2d, & |
---|
| 142 | ierr ) |
---|
[667] | 143 | |
---|
| 144 | ALLOCATE( total_2d_a(-nbgp:ny_a+nbgp,-nbgp:nx_a+nbgp), & |
---|
| 145 | total_2d_o(-nbgp:ny_o+nbgp,-nbgp:nx_o+nbgp) ) |
---|
| 146 | |
---|
| 147 | ENDIF |
---|
| 148 | |
---|
[108] | 149 | IF ( terminate_coupled_remote > 0 ) THEN |
---|
[274] | 150 | WRITE( message_string, * ) 'remote model "', & |
---|
| 151 | TRIM( coupling_mode_remote ), & |
---|
| 152 | '" terminated', & |
---|
| 153 | '&with terminate_coupled_remote = ', & |
---|
| 154 | terminate_coupled_remote, & |
---|
| 155 | '&local model "', TRIM( coupling_mode ), & |
---|
| 156 | '" has', & |
---|
| 157 | '&terminate_coupled = ', & |
---|
[667] | 158 | terminate_coupled |
---|
[258] | 159 | CALL message( 'surface_coupler', 'PA0310', 1, 2, 0, 6, 0 ) |
---|
[108] | 160 | RETURN |
---|
| 161 | ENDIF |
---|
[667] | 162 | |
---|
[291] | 163 | |
---|
[108] | 164 | ! |
---|
| 165 | !-- Exchange the current simulated time between the models, |
---|
[667] | 166 | !-- currently just for total_2ding |
---|
[709] | 167 | IF ( coupling_topology == 0 ) THEN |
---|
| 168 | |
---|
| 169 | CALL MPI_SEND( time_since_reference_point, 1, MPI_REAL, target_id, 11, & |
---|
| 170 | comm_inter, ierr ) |
---|
| 171 | CALL MPI_RECV( time_since_reference_point_rem, 1, MPI_REAL, target_id, & |
---|
| 172 | 11, comm_inter, status, ierr ) |
---|
[667] | 173 | ELSE |
---|
[709] | 174 | |
---|
[667] | 175 | IF ( myid == 0 ) THEN |
---|
[709] | 176 | |
---|
| 177 | CALL MPI_SEND( time_since_reference_point, 1, MPI_REAL, target_id, & |
---|
| 178 | 11, comm_inter, ierr ) |
---|
| 179 | CALL MPI_RECV( time_since_reference_point_rem, 1, MPI_REAL, & |
---|
[667] | 180 | target_id, 11, comm_inter, status, ierr ) |
---|
[709] | 181 | |
---|
[667] | 182 | ENDIF |
---|
[709] | 183 | |
---|
| 184 | CALL MPI_BCAST( time_since_reference_point_rem, 1, MPI_REAL, 0, comm2d, & |
---|
| 185 | ierr ) |
---|
| 186 | |
---|
[667] | 187 | ENDIF |
---|
[102] | 188 | |
---|
| 189 | ! |
---|
| 190 | !-- Exchange the interface data |
---|
| 191 | IF ( coupling_mode == 'atmosphere_to_ocean' ) THEN |
---|
[667] | 192 | |
---|
| 193 | ! |
---|
[709] | 194 | !-- Horizontal grid size and number of processors is equal in ocean and |
---|
| 195 | !-- atmosphere |
---|
| 196 | IF ( coupling_topology == 0 ) THEN |
---|
[102] | 197 | |
---|
| 198 | ! |
---|
[709] | 199 | !-- Send heat flux at bottom surface to the ocean |
---|
| 200 | CALL MPI_SEND( shf(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 12, & |
---|
| 201 | comm_inter, ierr ) |
---|
[102] | 202 | ! |
---|
[709] | 203 | !-- Send humidity flux at bottom surface to the ocean |
---|
[667] | 204 | IF ( humidity ) THEN |
---|
[709] | 205 | CALL MPI_SEND( qsws(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 13, & |
---|
| 206 | comm_inter, ierr ) |
---|
[667] | 207 | ENDIF |
---|
| 208 | ! |
---|
[709] | 209 | !-- Receive temperature at the bottom surface from the ocean |
---|
| 210 | CALL MPI_RECV( pt(0,nysg,nxlg), 1, type_xy, target_id, 14, & |
---|
| 211 | comm_inter, status, ierr ) |
---|
[108] | 212 | ! |
---|
[709] | 213 | !-- Send the momentum flux (u) at bottom surface to the ocean |
---|
| 214 | CALL MPI_SEND( usws(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 15, & |
---|
| 215 | comm_inter, ierr ) |
---|
[102] | 216 | ! |
---|
[709] | 217 | !-- Send the momentum flux (v) at bottom surface to the ocean |
---|
| 218 | CALL MPI_SEND( vsws(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 16, & |
---|
| 219 | comm_inter, ierr ) |
---|
[102] | 220 | ! |
---|
[709] | 221 | !-- Receive u at the bottom surface from the ocean |
---|
| 222 | CALL MPI_RECV( u(0,nysg,nxlg), 1, type_xy, target_id, 17, & |
---|
| 223 | comm_inter, status, ierr ) |
---|
[667] | 224 | ! |
---|
[709] | 225 | !-- Receive v at the bottom surface from the ocean |
---|
| 226 | CALL MPI_RECV( v(0,nysg,nxlg), 1, type_xy, target_id, 18, & |
---|
| 227 | comm_inter, status, ierr ) |
---|
[667] | 228 | ! |
---|
| 229 | !-- Horizontal grid size or number of processors differs between |
---|
| 230 | !-- ocean and atmosphere |
---|
| 231 | ELSE |
---|
| 232 | |
---|
| 233 | ! |
---|
[709] | 234 | !-- Send heat flux at bottom surface to the ocean |
---|
[1353] | 235 | total_2d_a = 0.0_wp |
---|
| 236 | total_2d = 0.0_wp |
---|
[667] | 237 | total_2d(nys:nyn,nxl:nxr) = shf(nys:nyn,nxl:nxr) |
---|
[709] | 238 | |
---|
| 239 | CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, MPI_SUM, 0, & |
---|
| 240 | comm2d, ierr ) |
---|
| 241 | CALL interpolate_to_ocean( 12 ) |
---|
[667] | 242 | ! |
---|
[709] | 243 | !-- Send humidity flux at bottom surface to the ocean |
---|
| 244 | IF ( humidity ) THEN |
---|
[1353] | 245 | total_2d_a = 0.0_wp |
---|
| 246 | total_2d = 0.0_wp |
---|
[667] | 247 | total_2d(nys:nyn,nxl:nxr) = qsws(nys:nyn,nxl:nxr) |
---|
[709] | 248 | |
---|
| 249 | CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, MPI_SUM, & |
---|
| 250 | 0, comm2d, ierr ) |
---|
| 251 | CALL interpolate_to_ocean( 13 ) |
---|
[667] | 252 | ENDIF |
---|
| 253 | ! |
---|
[709] | 254 | !-- Receive temperature at the bottom surface from the ocean |
---|
| 255 | IF ( myid == 0 ) THEN |
---|
[667] | 256 | CALL MPI_RECV( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, & |
---|
| 257 | target_id, 14, comm_inter, status, ierr ) |
---|
| 258 | ENDIF |
---|
| 259 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[709] | 260 | CALL MPI_BCAST( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, 0, comm2d, & |
---|
| 261 | ierr ) |
---|
[667] | 262 | pt(0,nysg:nyng,nxlg:nxrg) = total_2d_a(nysg:nyng,nxlg:nxrg) |
---|
| 263 | ! |
---|
[709] | 264 | !-- Send momentum flux (u) at bottom surface to the ocean |
---|
[1353] | 265 | total_2d_a = 0.0_wp |
---|
| 266 | total_2d = 0.0_wp |
---|
[667] | 267 | total_2d(nys:nyn,nxl:nxr) = usws(nys:nyn,nxl:nxr) |
---|
[709] | 268 | CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, MPI_SUM, 0, & |
---|
| 269 | comm2d, ierr ) |
---|
| 270 | CALL interpolate_to_ocean( 15 ) |
---|
[667] | 271 | ! |
---|
[709] | 272 | !-- Send momentum flux (v) at bottom surface to the ocean |
---|
[1353] | 273 | total_2d_a = 0.0_wp |
---|
| 274 | total_2d = 0.0_wp |
---|
[667] | 275 | total_2d(nys:nyn,nxl:nxr) = vsws(nys:nyn,nxl:nxr) |
---|
[709] | 276 | CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, MPI_SUM, 0, & |
---|
| 277 | comm2d, ierr ) |
---|
| 278 | CALL interpolate_to_ocean( 16 ) |
---|
[667] | 279 | ! |
---|
[709] | 280 | !-- Receive u at the bottom surface from the ocean |
---|
| 281 | IF ( myid == 0 ) THEN |
---|
[667] | 282 | CALL MPI_RECV( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, & |
---|
[709] | 283 | target_id, 17, comm_inter, status, ierr ) |
---|
[667] | 284 | ENDIF |
---|
| 285 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[709] | 286 | CALL MPI_BCAST( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, 0, comm2d, & |
---|
| 287 | ierr ) |
---|
[667] | 288 | u(0,nysg:nyng,nxlg:nxrg) = total_2d_a(nysg:nyng,nxlg:nxrg) |
---|
| 289 | ! |
---|
[709] | 290 | !-- Receive v at the bottom surface from the ocean |
---|
| 291 | IF ( myid == 0 ) THEN |
---|
[667] | 292 | CALL MPI_RECV( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, & |
---|
[709] | 293 | target_id, 18, comm_inter, status, ierr ) |
---|
[667] | 294 | ENDIF |
---|
| 295 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[709] | 296 | CALL MPI_BCAST( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, 0, comm2d, & |
---|
| 297 | ierr ) |
---|
[667] | 298 | v(0,nysg:nyng,nxlg:nxrg) = total_2d_a(nysg:nyng,nxlg:nxrg) |
---|
| 299 | |
---|
| 300 | ENDIF |
---|
| 301 | |
---|
[102] | 302 | ELSEIF ( coupling_mode == 'ocean_to_atmosphere' ) THEN |
---|
| 303 | |
---|
| 304 | ! |
---|
[667] | 305 | !-- Horizontal grid size and number of processors is equal |
---|
| 306 | !-- in ocean and atmosphere |
---|
| 307 | IF ( coupling_topology == 0 ) THEN |
---|
| 308 | ! |
---|
[709] | 309 | !-- Receive heat flux at the sea surface (top) from the atmosphere |
---|
| 310 | CALL MPI_RECV( tswst(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 12, & |
---|
| 311 | comm_inter, status, ierr ) |
---|
[102] | 312 | ! |
---|
[709] | 313 | !-- Receive humidity flux from the atmosphere (bottom) |
---|
[667] | 314 | !-- and add it to the heat flux at the sea surface (top)... |
---|
| 315 | IF ( humidity_remote ) THEN |
---|
| 316 | CALL MPI_RECV( qswst_remote(nysg,nxlg), ngp_xy, MPI_REAL, & |
---|
| 317 | target_id, 13, comm_inter, status, ierr ) |
---|
| 318 | ENDIF |
---|
| 319 | ! |
---|
| 320 | !-- Send sea surface temperature to the atmosphere model |
---|
[709] | 321 | CALL MPI_SEND( pt(nzt,nysg,nxlg), 1, type_xy, target_id, 14, & |
---|
| 322 | comm_inter, ierr ) |
---|
[667] | 323 | ! |
---|
| 324 | !-- Receive momentum flux (u) at the sea surface (top) from the atmosphere |
---|
[709] | 325 | CALL MPI_RECV( uswst(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 15, & |
---|
| 326 | comm_inter, status, ierr ) |
---|
[667] | 327 | ! |
---|
| 328 | !-- Receive momentum flux (v) at the sea surface (top) from the atmosphere |
---|
[709] | 329 | CALL MPI_RECV( vswst(nysg,nxlg), ngp_xy, MPI_REAL, target_id, 16, & |
---|
| 330 | comm_inter, status, ierr ) |
---|
[667] | 331 | ! |
---|
[709] | 332 | !-- Send u to the atmosphere |
---|
| 333 | CALL MPI_SEND( u(nzt,nysg,nxlg), 1, type_xy, target_id, 17, & |
---|
| 334 | comm_inter, ierr ) |
---|
[667] | 335 | ! |
---|
[709] | 336 | !-- Send v to the atmosphere |
---|
| 337 | CALL MPI_SEND( v(nzt,nysg,nxlg), 1, type_xy, target_id, 18, & |
---|
| 338 | comm_inter, ierr ) |
---|
| 339 | ! |
---|
[667] | 340 | !-- Horizontal gridsize or number of processors differs between |
---|
| 341 | !-- ocean and atmosphere |
---|
| 342 | ELSE |
---|
| 343 | ! |
---|
[709] | 344 | !-- Receive heat flux at the sea surface (top) from the atmosphere |
---|
| 345 | IF ( myid == 0 ) THEN |
---|
[667] | 346 | CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, & |
---|
[709] | 347 | target_id, 12, comm_inter, status, ierr ) |
---|
[667] | 348 | ENDIF |
---|
| 349 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[709] | 350 | CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, 0, comm2d, & |
---|
| 351 | ierr ) |
---|
[667] | 352 | tswst(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg) |
---|
| 353 | ! |
---|
[709] | 354 | !-- Receive humidity flux at the sea surface (top) from the atmosphere |
---|
| 355 | IF ( humidity_remote ) THEN |
---|
| 356 | IF ( myid == 0 ) THEN |
---|
[667] | 357 | CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, & |
---|
[709] | 358 | target_id, 13, comm_inter, status, ierr ) |
---|
[667] | 359 | ENDIF |
---|
| 360 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[709] | 361 | CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, 0, & |
---|
| 362 | comm2d, ierr) |
---|
[667] | 363 | qswst_remote(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg) |
---|
| 364 | ENDIF |
---|
| 365 | ! |
---|
| 366 | !-- Send surface temperature to atmosphere |
---|
[1353] | 367 | total_2d_o = 0.0_wp |
---|
| 368 | total_2d = 0.0_wp |
---|
[667] | 369 | total_2d(nys:nyn,nxl:nxr) = pt(nzt,nys:nyn,nxl:nxr) |
---|
| 370 | |
---|
[709] | 371 | CALL MPI_REDUCE( total_2d, total_2d_o, ngp_o, MPI_REAL, MPI_SUM, 0, & |
---|
| 372 | comm2d, ierr) |
---|
| 373 | CALL interpolate_to_atmos( 14 ) |
---|
[667] | 374 | ! |
---|
[709] | 375 | !-- Receive momentum flux (u) at the sea surface (top) from the atmosphere |
---|
| 376 | IF ( myid == 0 ) THEN |
---|
[667] | 377 | CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, & |
---|
[709] | 378 | target_id, 15, comm_inter, status, ierr ) |
---|
[667] | 379 | ENDIF |
---|
| 380 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 381 | CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, & |
---|
[709] | 382 | 0, comm2d, ierr ) |
---|
[667] | 383 | uswst(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg) |
---|
| 384 | ! |
---|
[709] | 385 | !-- Receive momentum flux (v) at the sea surface (top) from the atmosphere |
---|
| 386 | IF ( myid == 0 ) THEN |
---|
[667] | 387 | CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, & |
---|
[709] | 388 | target_id, 16, comm_inter, status, ierr ) |
---|
[667] | 389 | ENDIF |
---|
| 390 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[709] | 391 | CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, 0, comm2d, & |
---|
| 392 | ierr ) |
---|
[667] | 393 | vswst(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg) |
---|
| 394 | ! |
---|
| 395 | !-- Send u to atmosphere |
---|
[1353] | 396 | total_2d_o = 0.0_wp |
---|
| 397 | total_2d = 0.0_wp |
---|
[667] | 398 | total_2d(nys:nyn,nxl:nxr) = u(nzt,nys:nyn,nxl:nxr) |
---|
[709] | 399 | CALL MPI_REDUCE( total_2d, total_2d_o, ngp_o, MPI_REAL, MPI_SUM, 0, & |
---|
| 400 | comm2d, ierr ) |
---|
| 401 | CALL interpolate_to_atmos( 17 ) |
---|
[667] | 402 | ! |
---|
| 403 | !-- Send v to atmosphere |
---|
[1353] | 404 | total_2d_o = 0.0_wp |
---|
| 405 | total_2d = 0.0_wp |
---|
[667] | 406 | total_2d(nys:nyn,nxl:nxr) = v(nzt,nys:nyn,nxl:nxr) |
---|
[709] | 407 | CALL MPI_REDUCE( total_2d, total_2d_o, ngp_o, MPI_REAL, MPI_SUM, 0, & |
---|
| 408 | comm2d, ierr ) |
---|
| 409 | CALL interpolate_to_atmos( 18 ) |
---|
[667] | 410 | |
---|
| 411 | ENDIF |
---|
| 412 | |
---|
| 413 | ! |
---|
| 414 | !-- Conversions of fluxes received from atmosphere |
---|
| 415 | IF ( humidity_remote ) THEN |
---|
[108] | 416 | ! |
---|
[709] | 417 | !-- Here tswst is still the sum of atmospheric bottom heat fluxes, |
---|
| 418 | !-- * latent heat of vaporization in m2/s2, or 540 cal/g, or 40.65 kJ/mol |
---|
| 419 | !-- /(rho_atm(=1.0)*c_p) |
---|
[1427] | 420 | tswst = tswst + qswst_remote * l_v / cp |
---|
[709] | 421 | ! |
---|
[667] | 422 | !-- ...and convert it to a salinity flux at the sea surface (top) |
---|
[108] | 423 | !-- following Steinhorn (1991), JPO 21, pp. 1681-1683: |
---|
| 424 | !-- S'w' = -S * evaporation / ( rho_water * ( 1 - S ) ) |
---|
[1427] | 425 | saswst = -1.0_wp * sa(nzt,:,:) * 0.001_wp * qswst_remote / & |
---|
[2031] | 426 | ( rho_ocean(nzt,:,:) * ( 1.0_wp - sa(nzt,:,:) * 0.001_wp ) ) |
---|
[108] | 427 | ENDIF |
---|
| 428 | |
---|
| 429 | ! |
---|
[102] | 430 | !-- Adjust the kinematic heat flux with respect to ocean density |
---|
| 431 | !-- (constants are the specific heat capacities for air and water) |
---|
[667] | 432 | !-- now tswst is the ocean top heat flux |
---|
[2031] | 433 | tswst = tswst / rho_ocean(nzt,:,:) * cp / cpw |
---|
[102] | 434 | |
---|
| 435 | ! |
---|
[667] | 436 | !-- Adjust the momentum fluxes with respect to ocean density |
---|
[2031] | 437 | uswst = uswst / rho_ocean(nzt,:,:) |
---|
| 438 | vswst = vswst / rho_ocean(nzt,:,:) |
---|
[102] | 439 | |
---|
[667] | 440 | ENDIF |
---|
| 441 | |
---|
[709] | 442 | IF ( coupling_topology == 1 ) THEN |
---|
[667] | 443 | DEALLOCATE( total_2d_o, total_2d_a ) |
---|
| 444 | ENDIF |
---|
| 445 | |
---|
| 446 | CALL cpu_log( log_point(39), 'surface_coupler', 'stop' ) |
---|
| 447 | |
---|
| 448 | #endif |
---|
| 449 | |
---|
| 450 | END SUBROUTINE surface_coupler |
---|
| 451 | |
---|
| 452 | |
---|
| 453 | |
---|
[1682] | 454 | !------------------------------------------------------------------------------! |
---|
| 455 | ! Description: |
---|
| 456 | ! ------------ |
---|
| 457 | !> @todo Missing subroutine description. |
---|
| 458 | !------------------------------------------------------------------------------! |
---|
[709] | 459 | SUBROUTINE interpolate_to_atmos( tag ) |
---|
[667] | 460 | |
---|
[880] | 461 | #if defined( __parallel ) |
---|
| 462 | |
---|
[1320] | 463 | USE arrays_3d, & |
---|
| 464 | ONLY: total_2d_a, total_2d_o |
---|
[667] | 465 | |
---|
[1320] | 466 | USE indices, & |
---|
| 467 | ONLY: nbgp, nx, nx_a, nx_o, ny, ny_a, ny_o |
---|
| 468 | |
---|
| 469 | USE kinds |
---|
| 470 | |
---|
[1324] | 471 | USE pegrid |
---|
[1320] | 472 | |
---|
[667] | 473 | IMPLICIT NONE |
---|
| 474 | |
---|
[1682] | 475 | INTEGER(iwp) :: dnx !< |
---|
| 476 | INTEGER(iwp) :: dnx2 !< |
---|
| 477 | INTEGER(iwp) :: dny !< |
---|
| 478 | INTEGER(iwp) :: dny2 !< |
---|
| 479 | INTEGER(iwp) :: i !< |
---|
| 480 | INTEGER(iwp) :: ii !< |
---|
| 481 | INTEGER(iwp) :: j !< |
---|
| 482 | INTEGER(iwp) :: jj !< |
---|
[667] | 483 | |
---|
[1682] | 484 | INTEGER(iwp), intent(in) :: tag !< |
---|
[1320] | 485 | |
---|
[667] | 486 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 487 | |
---|
[709] | 488 | IF ( myid == 0 ) THEN |
---|
| 489 | ! |
---|
| 490 | !-- Cyclic boundary conditions for the total 2D-grid |
---|
[667] | 491 | total_2d_o(-nbgp:-1,:) = total_2d_o(ny+1-nbgp:ny,:) |
---|
| 492 | total_2d_o(:,-nbgp:-1) = total_2d_o(:,nx+1-nbgp:nx) |
---|
| 493 | |
---|
| 494 | total_2d_o(ny+1:ny+nbgp,:) = total_2d_o(0:nbgp-1,:) |
---|
| 495 | total_2d_o(:,nx+1:nx+nbgp) = total_2d_o(:,0:nbgp-1) |
---|
| 496 | |
---|
[102] | 497 | ! |
---|
[667] | 498 | !-- Number of gridpoints of the fine grid within one mesh of the coarse grid |
---|
| 499 | dnx = (nx_o+1) / (nx_a+1) |
---|
| 500 | dny = (ny_o+1) / (ny_a+1) |
---|
[102] | 501 | |
---|
| 502 | ! |
---|
[709] | 503 | !-- Distance for interpolation around coarse grid points within the fine |
---|
| 504 | !-- grid (note: 2*dnx2 must not be equal with dnx) |
---|
[667] | 505 | dnx2 = 2 * ( dnx / 2 ) |
---|
| 506 | dny2 = 2 * ( dny / 2 ) |
---|
[102] | 507 | |
---|
[1353] | 508 | total_2d_a = 0.0_wp |
---|
[102] | 509 | ! |
---|
[667] | 510 | !-- Interpolation from ocean-grid-layer to atmosphere-grid-layer |
---|
| 511 | DO j = 0, ny_a |
---|
| 512 | DO i = 0, nx_a |
---|
| 513 | DO jj = 0, dny2 |
---|
| 514 | DO ii = 0, dnx2 |
---|
| 515 | total_2d_a(j,i) = total_2d_a(j,i) & |
---|
| 516 | + total_2d_o(j*dny+jj,i*dnx+ii) |
---|
| 517 | ENDDO |
---|
| 518 | ENDDO |
---|
| 519 | total_2d_a(j,i) = total_2d_a(j,i) / ( ( dnx2 + 1 ) * ( dny2 + 1 ) ) |
---|
| 520 | ENDDO |
---|
| 521 | ENDDO |
---|
| 522 | ! |
---|
[709] | 523 | !-- Cyclic boundary conditions for atmosphere grid |
---|
[667] | 524 | total_2d_a(-nbgp:-1,:) = total_2d_a(ny_a+1-nbgp:ny_a,:) |
---|
| 525 | total_2d_a(:,-nbgp:-1) = total_2d_a(:,nx_a+1-nbgp:nx_a) |
---|
| 526 | |
---|
| 527 | total_2d_a(ny_a+1:ny_a+nbgp,:) = total_2d_a(0:nbgp-1,:) |
---|
| 528 | total_2d_a(:,nx_a+1:nx_a+nbgp) = total_2d_a(:,0:nbgp-1) |
---|
| 529 | ! |
---|
| 530 | !-- Transfer of the atmosphere-grid-layer to the atmosphere |
---|
[709] | 531 | CALL MPI_SEND( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, target_id, & |
---|
| 532 | tag, comm_inter, ierr ) |
---|
[102] | 533 | |
---|
| 534 | ENDIF |
---|
| 535 | |
---|
[667] | 536 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
[102] | 537 | |
---|
[880] | 538 | #endif |
---|
| 539 | |
---|
[667] | 540 | END SUBROUTINE interpolate_to_atmos |
---|
[102] | 541 | |
---|
[667] | 542 | |
---|
[1682] | 543 | !------------------------------------------------------------------------------! |
---|
| 544 | ! Description: |
---|
| 545 | ! ------------ |
---|
| 546 | !> @todo Missing subroutine description. |
---|
| 547 | !------------------------------------------------------------------------------! |
---|
[709] | 548 | SUBROUTINE interpolate_to_ocean( tag ) |
---|
[667] | 549 | |
---|
[880] | 550 | #if defined( __parallel ) |
---|
| 551 | |
---|
[1320] | 552 | USE arrays_3d, & |
---|
| 553 | ONLY: total_2d_a, total_2d_o |
---|
[667] | 554 | |
---|
[1320] | 555 | USE indices, & |
---|
| 556 | ONLY: nbgp, nx, nx_a, nx_o, ny, ny_a, ny_o |
---|
| 557 | |
---|
| 558 | USE kinds |
---|
| 559 | |
---|
[1324] | 560 | USE pegrid |
---|
[1320] | 561 | |
---|
[667] | 562 | IMPLICIT NONE |
---|
| 563 | |
---|
[1682] | 564 | INTEGER(iwp) :: dnx !< |
---|
| 565 | INTEGER(iwp) :: dny !< |
---|
| 566 | INTEGER(iwp) :: i !< |
---|
| 567 | INTEGER(iwp) :: ii !< |
---|
| 568 | INTEGER(iwp) :: j !< |
---|
| 569 | INTEGER(iwp) :: jj !< |
---|
| 570 | INTEGER(iwp), intent(in) :: tag !< |
---|
[667] | 571 | |
---|
[1682] | 572 | REAL(wp) :: fl !< |
---|
| 573 | REAL(wp) :: fr !< |
---|
| 574 | REAL(wp) :: myl !< |
---|
| 575 | REAL(wp) :: myr !< |
---|
[709] | 576 | |
---|
[667] | 577 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 578 | |
---|
[709] | 579 | IF ( myid == 0 ) THEN |
---|
[667] | 580 | |
---|
| 581 | ! |
---|
[709] | 582 | !-- Number of gridpoints of the fine grid within one mesh of the coarse grid |
---|
[667] | 583 | dnx = ( nx_o + 1 ) / ( nx_a + 1 ) |
---|
| 584 | dny = ( ny_o + 1 ) / ( ny_a + 1 ) |
---|
| 585 | |
---|
| 586 | ! |
---|
[709] | 587 | !-- Cyclic boundary conditions for atmosphere grid |
---|
[667] | 588 | total_2d_a(-nbgp:-1,:) = total_2d_a(ny+1-nbgp:ny,:) |
---|
| 589 | total_2d_a(:,-nbgp:-1) = total_2d_a(:,nx+1-nbgp:nx) |
---|
| 590 | |
---|
| 591 | total_2d_a(ny+1:ny+nbgp,:) = total_2d_a(0:nbgp-1,:) |
---|
| 592 | total_2d_a(:,nx+1:nx+nbgp) = total_2d_a(:,0:nbgp-1) |
---|
| 593 | ! |
---|
[709] | 594 | !-- Bilinear Interpolation from atmosphere grid-layer to ocean grid-layer |
---|
[667] | 595 | DO j = 0, ny |
---|
| 596 | DO i = 0, nx |
---|
| 597 | myl = ( total_2d_a(j+1,i) - total_2d_a(j,i) ) / dny |
---|
| 598 | myr = ( total_2d_a(j+1,i+1) - total_2d_a(j,i+1) ) / dny |
---|
| 599 | DO jj = 0, dny-1 |
---|
[709] | 600 | fl = myl*jj + total_2d_a(j,i) |
---|
| 601 | fr = myr*jj + total_2d_a(j,i+1) |
---|
[667] | 602 | DO ii = 0, dnx-1 |
---|
| 603 | total_2d_o(j*dny+jj,i*dnx+ii) = ( fr - fl ) / dnx * ii + fl |
---|
| 604 | ENDDO |
---|
| 605 | ENDDO |
---|
| 606 | ENDDO |
---|
| 607 | ENDDO |
---|
| 608 | ! |
---|
[709] | 609 | !-- Cyclic boundary conditions for ocean grid |
---|
[667] | 610 | total_2d_o(-nbgp:-1,:) = total_2d_o(ny_o+1-nbgp:ny_o,:) |
---|
| 611 | total_2d_o(:,-nbgp:-1) = total_2d_o(:,nx_o+1-nbgp:nx_o) |
---|
| 612 | |
---|
| 613 | total_2d_o(ny_o+1:ny_o+nbgp,:) = total_2d_o(0:nbgp-1,:) |
---|
| 614 | total_2d_o(:,nx_o+1:nx_o+nbgp) = total_2d_o(:,0:nbgp-1) |
---|
| 615 | |
---|
| 616 | CALL MPI_SEND( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, & |
---|
| 617 | target_id, tag, comm_inter, ierr ) |
---|
| 618 | |
---|
| 619 | ENDIF |
---|
| 620 | |
---|
| 621 | CALL MPI_BARRIER( comm2d, ierr ) |
---|
| 622 | |
---|
[880] | 623 | #endif |
---|
| 624 | |
---|
[667] | 625 | END SUBROUTINE interpolate_to_ocean |
---|