source: palm/trunk/SOURCE/surface_coupler.f90 @ 683

Last change on this file since 683 was 668, checked in by suehring, 14 years ago

last commit documented

  • Property svn:keywords set to Id
File size: 19.1 KB
Line 
1 SUBROUTINE surface_coupler
2
3!------------------------------------------------------------------------------!
4! Current revisions:
5! -----------------
6!
7! Former revisions:
8! ------------------
9! $Id: surface_coupler.f90 668 2010-12-23 13:22:58Z raasch $
10!
11! 667 2010-12-23 12:06:00Z suehring/gryschka
12! additional case for nonequivalent processor and grid topopolgy in ocean and
13! atmosphere added (coupling_topology = 1)
14! Added exchange of u and v from Ocean to Atmosphere
15!
16! 291 2009-04-16 12:07:26Z raasch
17! Coupling with independent precursor runs.
18! Output of messages replaced by message handling routine.
19!
20! 206 2008-10-13 14:59:11Z raasch
21! Implementation of a MPI-1 Coupling: replaced myid with target_id,
22! deleted __mpi2 directives
23!
24! 109 2007-08-28 15:26:47Z letzel
25! Initial revision
26!
27! Description:
28! ------------
29! Data exchange at the interface between coupled models
30!------------------------------------------------------------------------------!
31
32    USE arrays_3d
33    USE control_parameters
34    USE cpulog
35    USE grid_variables
36    USE indices
37    USE interfaces
38    USE pegrid
39
40    IMPLICIT NONE
41
42    INTEGER ::  i, j, k
43
44    REAL    ::  time_since_reference_point_rem
45    REAL    ::  total_2d(-nbgp:ny+nbgp,-nbgp:nx+nbgp)
46
47#if defined( __parallel )
48
49    CALL cpu_log( log_point(39), 'surface_coupler', 'start' )
50
51
52
53!
54!-- In case of model termination initiated by the remote model
55!-- (terminate_coupled_remote > 0), initiate termination of the local model.
56!-- The rest of the coupler must then be skipped because it would cause an MPI
57!-- intercomminucation hang.
58!-- If necessary, the coupler will be called at the beginning of the next
59!-- restart run.
60
61    IF ( coupling_topology == 0 ) THEN
62       CALL MPI_SENDRECV( terminate_coupled,        1, MPI_INTEGER, target_id,  &
63                          0,                                                    &
64                          terminate_coupled_remote, 1, MPI_INTEGER, target_id,  &
65                          0, comm_inter, status, ierr )
66    ELSE
67       IF ( myid == 0) THEN
68          CALL MPI_SENDRECV( terminate_coupled,        1, MPI_INTEGER, &
69                             target_id, 0,                             &
70                             terminate_coupled_remote, 1, MPI_INTEGER, & 
71                             target_id, 0,                             &
72                             comm_inter, status, ierr )
73       ENDIF
74       CALL MPI_BCAST( terminate_coupled_remote, 1, MPI_INTEGER, 0, comm2d, ierr)
75
76       ALLOCATE( total_2d_a(-nbgp:ny_a+nbgp,-nbgp:nx_a+nbgp),       &
77                 total_2d_o(-nbgp:ny_o+nbgp,-nbgp:nx_o+nbgp) )
78
79    ENDIF
80
81    IF ( terminate_coupled_remote > 0 )  THEN
82       WRITE( message_string, * ) 'remote model "',                         &
83                                  TRIM( coupling_mode_remote ),             &
84                                  '" terminated',                           &
85                                  '&with terminate_coupled_remote = ',      &
86                                  terminate_coupled_remote,                 &
87                                  '&local model  "', TRIM( coupling_mode ), &
88                                  '" has',                                  &
89                                  '&terminate_coupled = ',                  &
90                                   terminate_coupled
91       CALL message( 'surface_coupler', 'PA0310', 1, 2, 0, 6, 0 )
92       RETURN
93    ENDIF
94 
95
96!
97!-- Exchange the current simulated time between the models,
98!-- currently just for total_2ding
99    IF ( coupling_topology == 0 ) THEN   
100       CALL MPI_SEND( time_since_reference_point, 1, MPI_REAL, &
101                      target_id, 11, comm_inter, ierr )
102       CALL MPI_RECV( time_since_reference_point_rem, 1, MPI_REAL, &
103                      target_id, 11, comm_inter, status, ierr )
104    ELSE
105       IF ( myid == 0 ) THEN
106          CALL MPI_SEND( time_since_reference_point, 1, MPI_REAL, &
107                         target_id, 11, comm_inter, ierr )
108          CALL MPI_RECV( time_since_reference_point_rem, 1, MPI_REAL, &
109                         target_id, 11, comm_inter, status, ierr )
110       ENDIF
111       CALL MPI_BCAST( time_since_reference_point_rem, 1, MPI_REAL, &
112                       0, comm2d, ierr )
113    ENDIF
114    WRITE ( 9, * ) 'simulated time: ', simulated_time
115    WRITE ( 9, * ) 'time since start of coupling: ', &
116                  time_since_reference_point, ' remote: ', &
117                  time_since_reference_point_rem
118   CALL local_flush( 9 )
119 
120
121!
122!-- Exchange the interface data
123    IF ( coupling_mode == 'atmosphere_to_ocean' )  THEN
124   
125!
126!--    Horizontal grid size and number of processors is equal
127!--    in ocean and atmosphere
128       IF ( coupling_topology == 0 ) THEN
129
130!
131!--       Send heat flux at bottom surface to the ocean model
132          CALL MPI_SEND( shf(nysg,nxlg), ngp_xy, MPI_REAL, &
133                         target_id, 12, comm_inter, ierr )
134
135!
136!--       Send humidity flux at bottom surface to the ocean model
137          IF ( humidity )  THEN
138             CALL MPI_SEND( qsws(nysg,nxlg), ngp_xy, MPI_REAL, &
139                            target_id, 13, comm_inter, ierr )
140          ENDIF
141
142!
143!--       Receive temperature at the bottom surface from the ocean model
144          WRITE ( 9, * )  '*** receive pt from ocean'
145          CALL local_flush( 9 )
146          CALL MPI_RECV( pt(0,nysg,nxlg), 1, type_xy, &
147                         target_id, 14, comm_inter, status, ierr )
148
149!
150!--       Send the momentum flux (u) at bottom surface to the ocean model
151          CALL MPI_SEND( usws(nysg,nxlg), ngp_xy, MPI_REAL, &
152                         target_id, 15, comm_inter, ierr )
153
154!
155!--       Send the momentum flux (v) at bottom surface to the ocean model
156          CALL MPI_SEND( vsws(nysg,nxlg), ngp_xy, MPI_REAL, &
157                         target_id, 16, comm_inter, ierr )
158
159!
160!--       Receive u at the bottom surface from the ocean model
161          CALL MPI_RECV( u(0,nysg,nxlg), 1, type_xy, &
162                         target_id, 17, comm_inter, status, ierr )
163
164!
165!--       Receive v at the bottom surface from the ocean model
166          CALL MPI_RECV( v(0,nysg,nxlg), 1, type_xy, &
167                         target_id, 18,  comm_inter, status, ierr )
168
169!
170!--    Horizontal grid size or number of processors differs between
171!--    ocean and atmosphere
172       ELSE
173     
174!
175!--       Send heat flux at bottom surface to the ocean model
176          total_2d_a = 0.0
177          total_2d = 0.0
178          total_2d(nys:nyn,nxl:nxr) = shf(nys:nyn,nxl:nxr)
179          CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, &
180                           MPI_SUM, 0, comm2d, ierr )
181          CALL interpolate_to_ocean(12)
182   
183!
184!--       Send humidity flux at bottom surface to the ocean model
185          IF ( humidity ) THEN
186             total_2d_a = 0.0
187             total_2d = 0.0
188             total_2d(nys:nyn,nxl:nxr) = qsws(nys:nyn,nxl:nxr)
189             CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, &
190                              MPI_SUM, 0, comm2d, ierr )
191             CALL interpolate_to_ocean(13)
192          ENDIF
193
194!
195!--       Receive temperature at the bottom surface from the ocean model
196          IF ( myid == 0 ) THEN
197             CALL MPI_RECV( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
198                            target_id, 14, comm_inter, status, ierr )   
199          ENDIF
200          CALL MPI_BARRIER( comm2d, ierr )
201          CALL MPI_BCAST( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
202                          0, comm2d, ierr )
203          pt(0,nysg:nyng,nxlg:nxrg) = total_2d_a(nysg:nyng,nxlg:nxrg)
204
205!
206!--       Send momentum flux (u) at bottom surface to the ocean model
207          total_2d_a = 0.0 
208          total_2d = 0.0
209          total_2d(nys:nyn,nxl:nxr) = usws(nys:nyn,nxl:nxr)
210          CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, &
211                           MPI_SUM, 0, comm2d, ierr )
212          CALL interpolate_to_ocean(15)
213
214!
215!--       Send momentum flux (v) at bottom surface to the ocean model
216          total_2d_a = 0.0
217          total_2d = 0.0
218          total_2d(nys:nyn,nxl:nxr) = vsws(nys:nyn,nxl:nxr)
219          CALL MPI_REDUCE( total_2d, total_2d_a, ngp_a, MPI_REAL, &
220                           MPI_SUM, 0, comm2d, ierr )
221          CALL interpolate_to_ocean(16)
222
223!
224!--       Receive u at the bottom surface from the ocean model
225          IF ( myid == 0 ) THEN
226             CALL MPI_RECV( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
227                            target_id, 17, comm_inter, status, ierr )           
228          ENDIF
229          CALL MPI_BARRIER( comm2d, ierr )
230          CALL MPI_BCAST( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
231                          0, comm2d, ierr )
232          u(0,nysg:nyng,nxlg:nxrg) = total_2d_a(nysg:nyng,nxlg:nxrg)
233   
234!
235!--       Receive v at the bottom surface from the ocean model
236          IF ( myid == 0 ) THEN
237             CALL MPI_RECV( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
238                            target_id, 18, comm_inter, status, ierr )           
239          ENDIF
240          CALL MPI_BARRIER( comm2d, ierr )
241          CALL MPI_BCAST( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
242                          0, comm2d, ierr )
243          v(0,nysg:nyng,nxlg:nxrg) = total_2d_a(nysg:nyng,nxlg:nxrg)
244
245       ENDIF
246
247    ELSEIF ( coupling_mode == 'ocean_to_atmosphere' )  THEN
248
249!
250!--    Horizontal grid size and number of processors is equal
251!--    in ocean and atmosphere
252       IF ( coupling_topology == 0 ) THEN
253!
254!--       Receive heat flux at the sea surface (top) from the atmosphere model
255          CALL MPI_RECV( tswst(nysg,nxlg), ngp_xy, MPI_REAL, &
256                         target_id, 12, comm_inter, status, ierr )
257
258
259!
260!--       Receive humidity flux from the atmosphere model (bottom)
261!--       and add it to the heat flux at the sea surface (top)...
262          IF ( humidity_remote )  THEN
263             CALL MPI_RECV( qswst_remote(nysg,nxlg), ngp_xy, MPI_REAL, &
264                            target_id, 13, comm_inter, status, ierr )
265
266          ENDIF
267
268!
269!--       Send sea surface temperature to the atmosphere model
270          CALL MPI_SEND( pt(nzt,nysg,nxlg), 1, type_xy, &
271                         target_id, 14, comm_inter, ierr )
272
273!
274!--       Receive momentum flux (u) at the sea surface (top) from the atmosphere
275!--       model
276          WRITE ( 9, * )  '*** receive uswst from atmosphere'
277          CALL local_flush( 9 )
278          CALL MPI_RECV( uswst(nysg,nxlg), ngp_xy, MPI_REAL, &
279                         target_id, 15, comm_inter, status, ierr )
280
281!
282!--       Receive momentum flux (v) at the sea surface (top) from the atmosphere
283!--       model
284          CALL MPI_RECV( vswst(nysg,nxlg), ngp_xy, MPI_REAL, &
285                         target_id, 16, comm_inter, status, ierr )
286
287!--       Send u to the atmosphere model
288          CALL MPI_SEND( u(nzt,nysg,nxlg), 1, type_xy, &
289                         target_id, 17, comm_inter, ierr )
290
291!
292!--       Send v to the atmosphere model
293          CALL MPI_SEND( v(nzt,nysg,nxlg), 1, type_xy, &
294                         target_id, 18, comm_inter, ierr )
295
296!
297!--    Horizontal gridsize or number of processors differs between
298!--    ocean and atmosphere
299       ELSE
300
301!
302!--       Receive heat flux at the sea surface (top) from the atmosphere model
303          IF ( myid == 0 ) THEN
304             CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
305                            target_id, 12, comm_inter, status, ierr )           
306          ENDIF
307          CALL MPI_BARRIER( comm2d, ierr )
308          CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
309                          0, comm2d, ierr)
310          tswst(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg)
311
312!
313!--       Receive humidity flux at the sea surface (top) from the
314!--       atmosphere model
315          IF ( humidity_remote ) THEN
316             IF ( myid == 0 ) THEN
317                CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
318                               target_id, 13, comm_inter, status, ierr )           
319             ENDIF
320             CALL MPI_BARRIER( comm2d, ierr )
321             CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
322                             0, comm2d, ierr)
323             qswst_remote(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg)
324          ENDIF
325
326!
327!--       Send surface temperature to atmosphere
328          total_2d_o = 0.0
329          total_2d = 0.0
330          total_2d(nys:nyn,nxl:nxr) = pt(nzt,nys:nyn,nxl:nxr)
331
332          CALL MPI_REDUCE(total_2d, total_2d_o, ngp_o, &
333                          MPI_REAL, MPI_SUM, 0, comm2d, ierr) 
334
335          CALL interpolate_to_atmos(14)
336
337!
338!--       Receive momentum flux (u) at the sea surface (top) from the
339!--       atmosphere model
340          IF ( myid == 0 ) THEN
341             CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
342                            target_id, 15, comm_inter, status, ierr )           
343          ENDIF
344          CALL MPI_BARRIER( comm2d, ierr )
345          CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
346                          0, comm2d, ierr)
347          uswst(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg)
348
349!
350!--       Receive momentum flux (v) at the sea surface (top) from the
351!--       atmosphere model
352          IF ( myid == 0 ) THEN
353             CALL MPI_RECV( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
354                            target_id, 16, comm_inter, status, ierr )           
355          ENDIF
356          CALL MPI_BARRIER( comm2d, ierr )
357          CALL MPI_BCAST( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
358                          0, comm2d, ierr)
359          vswst(nysg:nyng,nxlg:nxrg) = total_2d_o(nysg:nyng,nxlg:nxrg)
360
361!
362!--       Send u to atmosphere
363          total_2d_o = 0.0 
364          total_2d = 0.0
365          total_2d(nys:nyn,nxl:nxr) = u(nzt,nys:nyn,nxl:nxr)
366          CALL MPI_REDUCE(total_2d, total_2d_o, ngp_o, MPI_REAL, &
367                          MPI_SUM, 0, comm2d, ierr) 
368          CALL interpolate_to_atmos(17)
369
370!
371!--       Send v to atmosphere
372          total_2d_o = 0.0
373          total_2d = 0.0
374          total_2d(nys:nyn,nxl:nxr) = v(nzt,nys:nyn,nxl:nxr)
375          CALL MPI_REDUCE(total_2d, total_2d_o, ngp_o, MPI_REAL, &
376                          MPI_SUM, 0, comm2d, ierr) 
377          CALL interpolate_to_atmos(18)
378       
379       ENDIF
380
381!
382!--    Conversions of fluxes received from atmosphere
383       IF ( humidity_remote )  THEN
384          !here tswst is still the sum of atmospheric bottom heat fluxes
385          tswst = tswst + qswst_remote * 2.2626108e6 / 1005.0
386          !*latent heat of vaporization in m2/s2, or 540 cal/g, or 40.65 kJ/mol
387          !/(rho_atm(=1.0)*c_p)
388!
389!--        ...and convert it to a salinity flux at the sea surface (top)
390!--       following Steinhorn (1991), JPO 21, pp. 1681-1683:
391!--       S'w' = -S * evaporation / ( rho_water * ( 1 - S ) )
392          saswst = -1.0 * sa(nzt,:,:) * qswst_remote /  &
393                    ( rho(nzt,:,:) * ( 1.0 - sa(nzt,:,:) ) )
394       ENDIF
395
396!
397!--    Adjust the kinematic heat flux with respect to ocean density
398!--    (constants are the specific heat capacities for air and water)
399!--    now tswst is the ocean top heat flux
400       tswst = tswst / rho(nzt,:,:) * 1005.0 / 4218.0
401
402!
403!--    Adjust the momentum fluxes with respect to ocean density
404       uswst = uswst / rho(nzt,:,:)
405       vswst = vswst / rho(nzt,:,:)
406
407
408    ENDIF
409
410    IF ( coupling_topology == 1 ) THEN
411       DEALLOCATE( total_2d_o, total_2d_a )
412    ENDIF
413
414    CALL cpu_log( log_point(39), 'surface_coupler', 'stop' )
415
416#endif
417
418  END SUBROUTINE surface_coupler
419
420
421
422  SUBROUTINE interpolate_to_atmos(tag)
423
424    USE arrays_3d
425    USE control_parameters
426    USE grid_variables
427    USE indices
428    USE pegrid
429
430    IMPLICIT NONE
431
432 
433    INTEGER             ::  dnx, dnx2, dny, dny2, i, ii, j, jj
434    INTEGER, intent(in) ::  tag
435
436    CALL MPI_BARRIER( comm2d, ierr )
437
438    IF ( myid == 0 ) THEN
439
440!
441!--    cyclic boundary conditions for the total 2D-grid
442       total_2d_o(-nbgp:-1,:) = total_2d_o(ny+1-nbgp:ny,:)
443       total_2d_o(:,-nbgp:-1) = total_2d_o(:,nx+1-nbgp:nx)
444
445       total_2d_o(ny+1:ny+nbgp,:) = total_2d_o(0:nbgp-1,:)
446       total_2d_o(:,nx+1:nx+nbgp) = total_2d_o(:,0:nbgp-1)
447
448!
449!--    Number of gridpoints of the fine grid within one mesh of the coarse grid
450       dnx = (nx_o+1) / (nx_a+1) 
451       dny = (ny_o+1) / (ny_a+1) 
452
453!
454!--    Distance for interpolation around coarse grid points within the fine grid
455!--    (note: 2*dnx2 must not be equal with dnx) 
456       dnx2 = 2 * ( dnx / 2 )
457       dny2 = 2 * ( dny / 2 )
458
459       total_2d_a = 0.0
460!
461!--    Interpolation from ocean-grid-layer to atmosphere-grid-layer
462       DO  j = 0, ny_a
463          DO  i = 0, nx_a 
464             DO  jj = 0, dny2
465                DO  ii = 0, dnx2
466                   total_2d_a(j,i) = total_2d_a(j,i) &
467                                     + total_2d_o(j*dny+jj,i*dnx+ii)
468                ENDDO
469             ENDDO
470             total_2d_a(j,i) = total_2d_a(j,i) / ( ( dnx2 + 1 ) * ( dny2 + 1 ) )
471          ENDDO
472       ENDDO
473!
474!--    cyclic boundary conditions for atmosphere grid
475       total_2d_a(-nbgp:-1,:) = total_2d_a(ny_a+1-nbgp:ny_a,:)
476       total_2d_a(:,-nbgp:-1) = total_2d_a(:,nx_a+1-nbgp:nx_a)
477       
478       total_2d_a(ny_a+1:ny_a+nbgp,:) = total_2d_a(0:nbgp-1,:)
479       total_2d_a(:,nx_a+1:nx_a+nbgp) = total_2d_a(:,0:nbgp-1)
480!
481!--    Transfer of the atmosphere-grid-layer to the atmosphere
482       CALL MPI_SEND( total_2d_a(-nbgp,-nbgp), ngp_a, MPI_REAL, &
483                      target_id, tag, comm_inter, ierr )
484
485    ENDIF
486
487    CALL MPI_BARRIER( comm2d, ierr )
488
489  END SUBROUTINE interpolate_to_atmos
490
491
492  SUBROUTINE interpolate_to_ocean(tag)
493
494    USE arrays_3d
495    USE control_parameters
496    USE grid_variables
497    USE indices
498    USE pegrid
499
500    IMPLICIT NONE
501
502    REAL                ::  fl, fr, myl, myr
503    INTEGER             ::  dnx, dny, i, ii, j, jj
504    INTEGER, intent(in) ::  tag
505
506    CALL MPI_BARRIER( comm2d, ierr )
507
508    IF ( myid == 0 ) THEN   
509
510!
511!      Number of gridpoints of the fine grid within one mesh of the coarse grid
512       dnx = ( nx_o + 1 ) / ( nx_a + 1 ) 
513       dny = ( ny_o + 1 ) / ( ny_a + 1 ) 
514
515!
516!--    cyclic boundary conditions for atmosphere grid
517       total_2d_a(-nbgp:-1,:) = total_2d_a(ny+1-nbgp:ny,:)
518       total_2d_a(:,-nbgp:-1) = total_2d_a(:,nx+1-nbgp:nx)
519       
520       total_2d_a(ny+1:ny+nbgp,:) = total_2d_a(0:nbgp-1,:)
521       total_2d_a(:,nx+1:nx+nbgp) = total_2d_a(:,0:nbgp-1)
522!
523!--    Bilinear Interpolation from atmosphere-grid-layer to ocean-grid-layer
524       DO  j = 0, ny
525          DO  i = 0, nx
526             myl = ( total_2d_a(j+1,i)   - total_2d_a(j,i)   ) / dny
527             myr = ( total_2d_a(j+1,i+1) - total_2d_a(j,i+1) ) / dny
528             DO  jj = 0, dny-1
529                fl = myl*jj  + total_2d_a(j,i) 
530                fr = myr*jj  + total_2d_a(j,i+1) 
531                DO  ii = 0, dnx-1
532                   total_2d_o(j*dny+jj,i*dnx+ii) = ( fr - fl ) / dnx * ii + fl
533                ENDDO
534             ENDDO
535          ENDDO
536       ENDDO
537!
538!--    cyclic boundary conditions for ocean grid
539       total_2d_o(-nbgp:-1,:) = total_2d_o(ny_o+1-nbgp:ny_o,:)
540       total_2d_o(:,-nbgp:-1) = total_2d_o(:,nx_o+1-nbgp:nx_o)
541
542       total_2d_o(ny_o+1:ny_o+nbgp,:) = total_2d_o(0:nbgp-1,:)
543       total_2d_o(:,nx_o+1:nx_o+nbgp) = total_2d_o(:,0:nbgp-1)
544       
545
546       CALL MPI_SEND( total_2d_o(-nbgp,-nbgp), ngp_o, MPI_REAL, &
547                      target_id, tag, comm_inter, ierr )
548
549    ENDIF
550
551    CALL MPI_BARRIER( comm2d, ierr ) 
552
553  END SUBROUTINE interpolate_to_ocean
Note: See TracBrowser for help on using the repository browser.