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

Last change on this file since 1462 was 1428, checked in by maronga, 10 years ago

last commit documented

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