WAVEWATCH III  beta 0.0.1
w3ogcmmd.F90
Go to the documentation of this file.
1 #include "w3macros.h"
2 !/ ------------------------------------------------------------------- /
3 MODULE w3ogcmmd
4  !/
5  !/ +-----------------------------------+
6  !/ | WAVEWATCH III NOAA/NCEP |
7  !/ | A. Thevenin |
8  !/ | FORTRAN 90 |
9  !/ | Last update : 22-Mar-2021 |
10  !/ +-----------------------------------+
11  !/
12  !/ Jul-2013 : Origination. ( version 4.18 )
13  !/ For upgrades see subroutines.
14  !/ Apr-2016 : Add comments (J. Pianezze) ( version 5.07 )
15  !/ 22-Mar-2021 : Add extra coupling variables ( version 7.13 )
16  !/
17  !/ Copyright 2009-2012 National Weather Service (NWS),
18  !/ National Oceanic and Atmospheric Administration. All rights
19  !/ reserved. WAVEWATCH III is a trademark of the NWS.
20  !/ No unauthorized use without permission.
21  !/
22  ! 1. Purpose :
23  !
24  ! Module used for coupling applications between oceanic model and WW3 with OASIS3-MCT
25  !
26  ! 2. Variables and types :
27  !
28  ! 3. Subroutines and functions :
29  !
30  ! Name Type Scope Description
31  ! ----------------------------------------------------------------
32  ! SND_FIELDS_TO_OCEAN Subr. Public Send fields to ocean model
33  ! RCV_FIELDS_FROM_OCEAN Subr. Public Receive fields from ocean model
34  ! ----------------------------------------------------------------
35  !
36  ! 4. Subroutines and functions used :
37  !
38  ! Name Type Module Description
39  ! ----------------------------------------------------------------
40  ! CPL_OASIS_SEND Subr. W3OACPMD Send fields
41  ! CPL_OASIS_RECV Subr. W3OACPMD Receive fields
42  ! ----------------------------------------------------------------
43  !
44  ! 5. Remarks
45  ! 6. Switches :
46  ! 7. Source code :
47  !
48  !/ ------------------------------------------------------------------- /
49  !
50  IMPLICIT NONE
51  !
52  include "mpif.h"
53  !
54  PRIVATE
55  !
56  ! * Accessibility
57  PUBLIC snd_fields_to_ocean
59  !
60 CONTAINS
61  !/ ------------------------------------------------------------------- /
62  SUBROUTINE snd_fields_to_ocean()
63  !/
64  !/ +-----------------------------------+
65  !/ | WAVEWATCH III NOAA/NCEP |
66  !/ | A. Thevenin |
67  !/ | FORTRAN 90 |
68  !/ | Last update : 22-Mar-2021 |
69  !/ +-----------------------------------+
70  !/
71  !/ Jul-2013 : Origination. ( version 4.18 )
72  !/ Apr-2016 : Add comments (J. Pianezze) ( version 5.07 )
73  !/ 22-Mar-2021 : Add extra coupling variables ( version 7.13 )
74  !/
75  ! 1. Purpose :
76  !
77  ! Send coupling fields to oceanic model
78  !
79  ! 2. Method :
80  ! 3. Parameters :
81  ! 4. Subroutines used :
82  !
83  ! Name Type Module Description
84  ! -------------------------------------------------------------------
85  ! CPL_OASIS_SND Subr. W3OACPMD Send field to atmos/ocean model
86  ! -------------------------------------------------------------------
87  !
88  ! 5. Called by :
89  !
90  ! Name Type Module Description
91  ! ------------------------------------------------------------------
92  ! W3WAVE Subr. W3WAVEMD Wave model
93  ! ------------------------------------------------------------------
94  !
95  ! 6. Error messages :
96  ! 7. Remarks :
97  !
98  ! According to the present implementation, fields are sent at each coupling time step to OASIS
99  ! Consequently, OASIS cannot estimate any time average
100  ! For such an application, one must estimate the fields at each time step
101  ! (or a time step smaller than the coupling time step).
102  ! In such conditions, OASIS get the information every time step
103  ! but only send the information to the other code when the time matches the coupling time.
104  !
105  ! 8. Structure :
106  ! 9. Switches :
107  ! 10. Source code :
108  !
109  !/ ------------------------------------------------------------------- /
110  !
112  USE w3gdatmd, ONLY: nseal, mapsta, mapsf
113  USE w3adatmd, ONLY: hs, t0m1, t01, thm, bhd, tauox, tauoy, phioc,&
114  uba, ubd, tauwix, tauwiy, tusx, tusy, ussx, &
116  tauocy, wnmean
117  USE w3odatmd, ONLY: naproc, iaproc, undef
118  USE constants, ONLY: pi, dera
119  !
120  !/ ------------------------------------------------------------------- /
121  !/ Parameter list
122  !/
123  !
124  !/ ------------------------------------------------------------------- /
125  !/ Local parameters
126  !/
127  INTEGER :: i, isea, ix, iy
128  INTEGER, DIMENSION(NSEAL) :: mask
129  REAL(kind=8), dimension(nseal,1) :: rla_oasis_snd
130  INTEGER :: ib_do
131  LOGICAL :: ll_action
132  REAL(kind=8), dimension(nseal) :: tmp
133  !
134  !----------------------------------------------------------------------
135  ! * Executable part
136  !
137  DO i = 1, nseal
138  isea = iaproc + (i-1)*naproc
139  ix = mapsf(isea,1)
140  iy = mapsf(isea,2)
141  ! Get the mask : 1 - sea 0 - open boundary cells dried cells
142  mask(i) = mod(mapsta(iy,ix),2)
143  END DO
144  !
145  DO ib_do = 1, il_nb_snd
146  !
147  ! Mask - wet-drying
148  ! ---------------------------------------------------------------------
149  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_ODRY') THEN
150  tmp(1:nseal) = 0.0
151  WHERE(mask(1:nseal) /= undef) tmp(1:nseal)=mask(1:nseal)
152  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
153  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
154  ENDIF
155  !
156  ! Mean wave period (tmn in s) (m0,-1)
157  ! ---------------------------------------------------------------------
158  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_T0M1') THEN
159  tmp(1:nseal) = 0.0
160  WHERE(t0m1(1:nseal) /= undef) tmp(1:nseal)=t0m1(1:nseal)
161  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
162  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
163  ENDIF
164  !
165  ! Mean wave period (tmn in s) (m0,1)
166  ! ---------------------------------------------------------------------
167  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__T01') THEN
168  tmp(1:nseal) = 0.0
169  WHERE(t01(1:nseal) /= undef) tmp(1:nseal)=t01(1:nseal)
170  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
171  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
172  ENDIF
173  !
174  ! Mean wave number (wnm in m-1)
175  ! ---------------------------------------------------------------------
176  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__WNM') THEN
177  tmp(1:nseal) = 0.0
178  WHERE(wnmean(1:nseal) /= undef) tmp(1:nseal)=wnmean(1:nseal)
179  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
180  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
181  ENDIF
182  !
183  ! Charnock coefficient (-)
184  ! ---------------------------------------------------------------------
185  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_OCHA') THEN
186  tmp(1:nseal) = 0.0
187  WHERE(charn(1:nseal) /= undef) tmp(1:nseal)=charn(1:nseal)
188  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
189  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
190  ENDIF
191  !
192  ! Wave height (hs in m)
193  ! ---------------------------------------------------------------------
194  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__OHS') THEN
195  tmp(1:nseal) = 0.0
196  WHERE(hs(1:nseal) /= undef) tmp(1:nseal)=hs(1:nseal)
197  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
198  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
199  ENDIF
200  !
201  ! Cosinus of Mean wave direction (cos(theta) in radians)
202  ! ---------------------------------------------------------------------
203  ! dir : nautical convention (GRIDDED files) - 0 degree from north, 90 from east
204  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_CDIR') THEN
205  tmp(1:nseal) = 0.0
206  WHERE(thm(1:nseal) /= undef) tmp(1:nseal)=cos(thm(1:nseal))
207  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
208  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
209  ENDIF
210  !
211  ! Sinus of Mean wave direction (sin(theta) in radians)
212  ! ---------------------------------------------------------------------
213  ! dir : nautical convention (GRIDDED files) - 0 degree from north, 90 from east
214  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_SDIR') THEN
215  tmp(1:nseal) = 0.0
216  WHERE(thm(1:nseal) /= undef) tmp(1:nseal)=sin(thm(1:nseal))
217  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
218  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
219  ENDIF
220  !
221  ! Mean wave direction theta in radians
222  ! ---------------------------------------------------------------------
223  ! dir : nautical convention (GRIDDED files) - 0 degree from north, 90 from east
224  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__DIR') THEN
225  tmp(1:nseal) = 0.0
226  WHERE(thm /= undef) tmp(1:nseal)=thm(1:nseal)
227  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
228  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
229  ENDIF
230  !
231  ! Wave-induced Bernoulli head pressure (bhd in N.m-1) (J term, Smith JPO 2006)
232  ! ---------------------------------------------------------------------
233  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__BHD') THEN
234  tmp(1:nseal) = 0.0
235  WHERE(bhd(1:nseal) /= undef) tmp(1:nseal)=bhd(1:nseal)
236  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
237  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
238  ENDIF
239  !
240  ! Wave-ocean momentum flux (tauox in m2.s-2)
241  ! ---------------------------------------------------------------------
242  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TWOX') THEN
243  tmp(1:nseal) = 0.0
244  WHERE(tauox(1:nseal) /= undef) tmp(1:nseal)=tauox(1:nseal)
245  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
246  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
247  ENDIF
248  !
249  ! Wave-ocean momentum flux (tauoy in m2.s-2)
250  ! ---------------------------------------------------------------------
251  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TWOY') THEN
252  tmp(1:nseal) = 0.0
253  WHERE(tauoy(1:nseal) /= undef) tmp(1:nseal)=tauoy(1:nseal)
254  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
255  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
256  ENDIF
257  !
258  ! Wave-ocean total momentum flux (tauocx in Pa)
259  ! ---------------------------------------------------------------------
260  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TOCX') THEN
261  tmp(1:nseal) = 0.0
262  WHERE(tauocx(1:nseal) /= undef) tmp(1:nseal)=tauocx(1:nseal)
263  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
264  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
265  ENDIF
266  !
267  ! Wave-ocean total momentum flux (tauocy in Pa)
268  ! ---------------------------------------------------------------------
269  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TOCY') THEN
270  tmp(1:nseal) = 0.0
271  WHERE(tauocy(1:nseal) /= undef) tmp(1:nseal)=tauocy(1:nseal)
272  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
273  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
274  ENDIF
275  !
276  ! Wave-to-ocean TKE flux (phioc in W.m-2)
277  ! ---------------------------------------------------------------------
278  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__FOC') THEN
279  tmp(1:nseal) = 0.0
280  WHERE(phioc(1:nseal) /= undef) tmp(1:nseal)=phioc(1:nseal)
281  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
282  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
283  ENDIF
284  !
285  ! Momentum flux due to bottom friction (taubblx in m2.s-2)
286  ! ---------------------------------------------------------------------
287  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TBBX') THEN
288  tmp(1:nseal) = 0.0
289  WHERE(taubbl(1:nseal,1) /= undef) tmp(1:nseal)=taubbl(1:nseal,1)
290  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
291  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
292  ENDIF
293  !
294  ! Momentum flux due to bottom friction (taubbly in m2.s-2)
295  ! ---------------------------------------------------------------------
296  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TBBY') THEN
297  tmp(1:nseal) = 0.0
298  WHERE(taubbl(1:nseal,2) /= undef) tmp(1:nseal)=taubbl(1:nseal,2)
299  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
300  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
301  ENDIF
302  !
303  ! Energy flux due to bottom friction (phibbl in W.m-2)
304  ! ---------------------------------------------------------------------
305  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__FBB') THEN
306  tmp(1:nseal) = 0.0
307  WHERE(phibbl(1:nseal) /= undef) tmp(1:nseal)=phibbl(1:nseal)
308  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
309  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
310  ENDIF
311  !
312  ! rms amplitude of orbital velocity of the waves (ubr in m.s-1)
313  ! ---------------------------------------------------------------------
314  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__UBR') THEN
315  tmp(1:nseal) = 0.0
316  WHERE(uba(1:nseal) /= undef) tmp(1:nseal)=uba(1:nseal)
317  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
318  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
319  ENDIF
320  !
321  ! x component of the near-bottom rms wave velocity (in m.s-1)
322  ! ---------------------------------------------------------------------
323  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_UBRX') THEN
324  tmp(1:nseal) = 0.0
325  WHERE(uba(1:nseal) /= undef) tmp(1:nseal)=uba(1:nseal)*cos(ubd(1:nseal))
326  rla_oasis_snd(:,1) = tmp(1:nseal)
327  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
328  ENDIF
329  !
330  ! y component of the near-bottom rms wave velocity (in m.s-1)
331  ! ---------------------------------------------------------------------
332  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_UBRY') THEN
333  tmp(1:nseal) = 0.0
334  WHERE(uba(1:nseal) /= undef) tmp(1:nseal)=uba(1:nseal)*sin(ubd(1:nseal))
335  rla_oasis_snd(:,1) = tmp(1:nseal)
336  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
337  ENDIF
338  !
339  ! Net wave-supported stress, u component (tauwix in m2.s-2)
340  ! ---------------------------------------------------------------------
341  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TAWX') THEN
342  tmp(1:nseal) = 0.0
343  WHERE(tauwix(1:nseal) /= undef) tmp(1:nseal)=tauwix(1:nseal)
344  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
345  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
346  ENDIF
347  !
348  ! Net wave-supported stress, v component (tauwix in m2.s-2)
349  ! ---------------------------------------------------------------------
350  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TAWY') THEN
351  tmp(1:nseal) = 0.0
352  WHERE(tauwiy(1:nseal) /= undef) tmp(1:nseal)=tauwiy(1:nseal)
353  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
354  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
355  ENDIF
356  !
357  ! Volume transport associated to Stokes drift, u component (tusx in m2.s-1)
358  ! ---------------------------------------------------------------------
359  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TUSX') THEN
360  tmp(1:nseal) = 0.0
361  WHERE(tusx(1:nseal) /= undef) tmp(1:nseal)=tusx(1:nseal)
362  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
363  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
364  ENDIF
365  !
366  ! Volume transport associated to Stokes drift, v component (tusy in m2.s-1)
367  ! ---------------------------------------------------------------------
368  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_TUSY') THEN
369  tmp(1:nseal) = 0.0
370  WHERE(tusy(1:nseal) /= undef) tmp(1:nseal)=tusy(1:nseal)
371  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
372  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
373  ENDIF
374  !
375  ! Surface Stokes drift, u component (ussx in m.s-1)
376  ! ---------------------------------------------------------------------
377  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_USSX') THEN
378  tmp(1:nseal) = 0.0
379  WHERE(ussx(1:nseal) /= undef) tmp(1:nseal)=ussx(1:nseal)
380  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
381  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
382  ENDIF
383  !
384  ! Surface Stokes drift, v component (ussy in m.s-1)
385  ! ---------------------------------------------------------------------
386  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_USSY') THEN
387  tmp(1:nseal) = 0.0
388  WHERE(ussy(1:nseal) /= undef) tmp(1:nseal)=ussy(1:nseal)
389  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
390  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
391  ENDIF
392  !
393  ! Mean wave length (wlm in m)
394  ! ---------------------------------------------------------------------
395  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3___LM') THEN
396  tmp(1:nseal) = 0.0
397  WHERE(wlm(1:nseal) /= undef) tmp(1:nseal)=wlm(1:nseal)
398  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
399  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
400  ENDIF
401  !
402  ENDDO
403  !/ ------------------------------------------------------------------- /
404  END SUBROUTINE snd_fields_to_ocean
405  !/ ------------------------------------------------------------------- /
406  SUBROUTINE rcv_fields_from_ocean(ID_LCOMM, IDFLD, FXN, FYN, FAN)
407  !/
408  !/ +-----------------------------------+
409  !/ | WAVEWATCH III NOAA/NCEP |
410  !/ | A. Thevenin |
411  !/ | FORTRAN 90 |
412  !/ | Last update : April-2016 |
413  !/ +-----------------------------------+
414  !/
415  !/ Jul-2013 : Origination. ( version 4.18 )
416  !/ Apr-2014 : Add IDFLD, FXN, FYX and FAN (M. Accensi) ( version 5.07 )
417  !/ Apr-2016 : Add comments (J. Pianezze) ( version 5.07 )
418  !/
419  ! 1. Purpose :
420  !
421  ! Receive coupling fields from oceanic model
422  !
423  ! 2. Method :
424  ! 3. Parameters :
425  !
426  ! Parameter list
427  ! ----------------------------------------------------------------
428  ! ID_LCOMM Char. I MPI communicator
429  ! IDFLD Int. I Name of the exchange fields
430  ! FXN Int. I/O First exchange field
431  ! FYN Int. I/O Second exchange field
432  ! FAN Int. I/O Third exchange field
433  ! ----------------------------------------------------------------
434  !
435  ! 4. Subroutines used :
436  !
437  ! Name Type Module Description
438  ! -------------------------------------------------------------------
439  ! CPL_OASIS_RCV Subr. W3OACPMD Receive fields from atmos/ocean model
440  ! W3S2XY Subr. W3SERVMD Convert from storage (NSEA) to spatial grid (NX, NY)
441  ! -------------------------------------------------------------------
442  !
443  ! 5. Called by :
444  !
445  ! Name Type Module Description
446  ! ------------------------------------------------------------------
447  ! W3FLDG Subr. W3FLDSMD Manage input fields of depth,
448  ! current, wind and ice concentration
449  ! ------------------------------------------------------------------
450  !
451  ! 6. Error messages :
452  ! 7. Remarks :
453  !
454  ! IDFLD C*3 I/O ID string for field type, valid are: 'LEV', 'CUR' (J=1,2)
455  !
456  ! 8. Structure :
457  ! 9. Switches :
458  ! 10. Source code :
459  !
460  !/ ------------------------------------------------------------------- /
461  !
463  USE w3gdatmd, ONLY: nx, ny, nseal, nsea, mapsf
464  USE w3odatmd, ONLY: naproc, iaproc
465  USE w3servmd, ONLY: w3s2xy
466  !
467  !/ ------------------------------------------------------------------- /
468  !/ Parameter list
469  !/
470  INTEGER, INTENT(IN) :: id_lcomm
471  CHARACTER(LEN=3), INTENT(IN) :: idfld
472  REAL, INTENT(INOUT) :: fxn(:,:), fyn(:,:), fan(:,:)
473  !
474  !/ ------------------------------------------------------------------- /
475  !/ Local parameters
476  !/
477  LOGICAL :: ll_action
478  INTEGER :: ib_do, ib_i, ib_j, il_err
479  INTEGER, SAVE :: id_oasis_time_wetdryonlyonce = -1
480  REAL(kind=8), dimension(nseal,1) :: rla_oasis_rcv
481  REAL(kind=8), dimension(nseal) :: tmp, maskt, masku, maskv
482  REAL, DIMENSION(1:NSEA) :: snd_buff,rcv_buff
483  !
484  !----------------------------------------------------------------------
485  ! * Executable part
486  !
487  maskt(:)=1.
488  masku(:)=1.
489  maskv(:)=1.
490  rla_oasis_rcv(:,:) = 0.0
491  !
492  ! ---------------------------------------------------------------------
493  ! Perform mask variables
494  ! ---------------------------------------------------------------------
495  !
496  ! For the same coupling time, W3FLDG is called for the level and current variables.
497  ! As RCV_FIELDS_FROM_OCEAN is called from W3FLDG, the following test prevents to
498  ! exchange the wet-dry variables more than once per coupling time.
499  !cval well but it cannot work because MASKT,MASKU,MASKV variable are not global variable
500  !cval Anyway we will give up the exchange of mask, it is not a good idea at all
501 
502  IF (id_oasis_time > id_oasis_time_wetdryonlyonce) THEN
503  !
504  DO ib_do = 1, il_nb_rcv
505  !
506  ! Land mask - u
507  ! ---------------------------------------------------------------------
508  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_OWDH') THEN
509 
510  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
511  IF (ll_action) THEN
512  maskt(1:nseal) = rla_oasis_rcv(1:nseal,1)
513  ENDIF
514  ENDIF
515  !
516  ! Land mask - h
517  ! ---------------------------------------------------------------------
518  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_OWDU') THEN
519  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
520  IF (ll_action) THEN
521  masku(1:nseal) = rla_oasis_rcv(1:nseal,1)
522  ENDIF
523  ENDIF
524  !
525  ! Land mask - v
526  ! ---------------------------------------------------------------------
527  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_OWDV') THEN
528  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
529  IF (ll_action) THEN
530  maskv(1:nseal) = rla_oasis_rcv(1:nseal,1)
531  ENDIF
532  ENDIF
533  !
534  ENDDO
535  !
536  ENDIF
537  !
538  ! ---------------------------------------------------------------------
539  ! Treatment of the dynamical variables
540  ! ---------------------------------------------------------------------
541  DO ib_do = 1, il_nb_rcv
542  !
543  ! Sea surface Height (m)
544  ! ---------------------------------------------------------------------
545  IF (idfld == 'LEV') THEN
546  !
547  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3__SSH') THEN
548  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
549  IF (ll_action) THEN
550  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1) * maskt(1:nseal)
551  snd_buff(1:nsea) = 0.0
552  DO ib_i = 1, nseal
553  ib_j = iaproc + (ib_i-1)*naproc
554  snd_buff(ib_j) = tmp(ib_i)
555  ENDDO
556  !
557  CALL mpi_allreduce(snd_buff(1:nsea), &
558  rcv_buff(1:nsea), &
559  nsea, &
560  mpi_real, &
561  mpi_sum, &
562  id_lcomm, &
563  il_err)
564  !
565  ! Convert from storage (NSEA) to spatial grid (NX, NY)
566  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fan)
567  !
568  ENDIF
569  ENDIF
570  ENDIF
571  !
572  ! Ocean sea surface current (m.s-1)
573  ! ---------------------------------------------------------------------
574  IF (idfld == 'CUR') THEN
575  !
576  ! u-component
577  ! ---------------------------------------------------------------------
578  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_OSSU') THEN
579  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
580  IF (ll_action) THEN
581  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1) * masku(1:nseal)
582  snd_buff(1:nsea) = 0.0
583  DO ib_i = 1, nseal
584  ib_j = iaproc + (ib_i-1)*naproc
585  snd_buff(ib_j) = tmp(ib_i)
586  ENDDO
587  !
588  CALL mpi_allreduce(snd_buff(1:nsea), &
589  rcv_buff(1:nsea), &
590  nsea, &
591  mpi_real, &
592  mpi_sum, &
593  id_lcomm, &
594  il_err)
595  !
596  ! Convert from storage (NSEA) to spatial grid (NX, NY)
597  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fxn)
598  !
599  ENDIF
600  ENDIF
601  !
602  ! v-component
603  ! ---------------------------------------------------------------------
604  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_OSSV') THEN
605  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
606  IF (ll_action) THEN
607  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1) * maskv(1:nseal)
608  snd_buff(1:nsea) = 0.0
609  DO ib_i = 1, nseal
610  ib_j = iaproc + (ib_i-1)*naproc
611  snd_buff(ib_j) = tmp(ib_i)
612  ENDDO
613  !
614  CALL mpi_allreduce(snd_buff(1:nsea), &
615  rcv_buff(1:nsea), &
616  nsea, &
617  mpi_real, &
618  mpi_sum, &
619  id_lcomm, &
620  il_err)
621  !
622  ! Convert from storage (NSEA) to spatial grid (NX, NY)
623  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fyn)
624  !
625  ENDIF
626  ENDIF
627  ENDIF
628  ENDDO
629  !
630  id_oasis_time_wetdryonlyonce = id_oasis_time
631  !
632  !/ ------------------------------------------------------------------- /
633  END SUBROUTINE rcv_fields_from_ocean
634  !/ ------------------------------------------------------------------- /
635  !/
636 END MODULE w3ogcmmd
637 !/
638 !/ ------------------------------------------------------------------- /
w3gdatmd::nseal
integer, pointer nseal
Definition: w3gdatmd.F90:1097
constants::pi
real, parameter pi
PI Value of Pi.
Definition: constants.F90:71
include
cmake src_list cmake include(${CMAKE_CURRENT_SOURCE_DIR}/cmake/check_switches.cmake) check_switches("$
Definition: CMakeLists.txt:15
w3adatmd::charn
real, dimension(:), pointer charn
Definition: w3adatmd.F90:603
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3adatmd::ussy
real, dimension(:), pointer ussy
Definition: w3adatmd.F90:607
constants::dera
real, parameter dera
DERA Conversion factor from degrees to radians.
Definition: constants.F90:77
w3adatmd::tusx
real, dimension(:), pointer tusx
Definition: w3adatmd.F90:607
w3oacpmd::cpl_oasis_snd
subroutine, public cpl_oasis_snd(ID_NB, ID_TIME, RDA_FIELD, LD_ACTION)
Definition: w3oacpmd.F90:563
w3adatmd::tusy
real, dimension(:), pointer tusy
Definition: w3adatmd.F90:607
w3oacpmd::snd_fld
type(cpl_field), dimension(ip_maxfld), public snd_fld
Definition: w3oacpmd.F90:76
w3oacpmd::il_nb_rcv
integer, public il_nb_rcv
Definition: w3oacpmd.F90:67
w3adatmd::t01
real, dimension(:), pointer t01
Definition: w3adatmd.F90:587
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3adatmd::tauocy
real, dimension(:), pointer tauocy
Definition: w3adatmd.F90:607
w3adatmd::hs
real, dimension(:), pointer hs
Definition: w3adatmd.F90:587
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
w3oacpmd::il_nb_snd
integer, public il_nb_snd
Definition: w3oacpmd.F90:67
w3adatmd::uba
real, dimension(:), pointer uba
Definition: w3adatmd.F90:614
w3adatmd::tauwix
real, dimension(:), pointer tauwix
Definition: w3adatmd.F90:603
w3adatmd::tauwiy
real, dimension(:), pointer tauwiy
Definition: w3adatmd.F90:603
w3adatmd::taubbl
real, dimension(:,:), pointer taubbl
Definition: w3adatmd.F90:614
w3adatmd::phibbl
real, dimension(:), pointer phibbl
Definition: w3adatmd.F90:614
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3ogcmmd
Definition: w3ogcmmd.F90:3
w3odatmd
Definition: w3odatmd.F90:3
w3oacpmd
Definition: w3oacpmd.F90:3
w3adatmd::bhd
real, dimension(:), pointer bhd
Definition: w3adatmd.F90:607
w3adatmd::wlm
real, dimension(:), pointer wlm
Definition: w3adatmd.F90:587
w3gdatmd::mapsf
integer, dimension(:,:), pointer mapsf
Definition: w3gdatmd.F90:1163
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
w3adatmd::wnmean
real, dimension(:), pointer wnmean
Definition: w3adatmd.F90:587
w3oacpmd::rcv_fld
type(cpl_field), dimension(ip_maxfld), public rcv_fld
Definition: w3oacpmd.F90:76
w3adatmd::phioc
real, dimension(:), pointer phioc
Definition: w3adatmd.F90:607
w3oacpmd::id_oasis_time
integer, public id_oasis_time
Definition: w3oacpmd.F90:78
w3oacpmd::cpl_oasis_rcv
subroutine, public cpl_oasis_rcv(ID_NB, ID_TIME, RDA_FIELD, LD_ACTION)
Definition: w3oacpmd.F90:634
w3ogcmmd::snd_fields_to_ocean
subroutine, public snd_fields_to_ocean()
Definition: w3ogcmmd.F90:63
w3adatmd::tauox
real, dimension(:), pointer tauox
Definition: w3adatmd.F90:607
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3adatmd::tauoy
real, dimension(:), pointer tauoy
Definition: w3adatmd.F90:607
w3gdatmd
Definition: w3gdatmd.F90:16
w3adatmd::ussx
real, dimension(:), pointer ussx
Definition: w3adatmd.F90:607
w3ogcmmd::rcv_fields_from_ocean
subroutine, public rcv_fields_from_ocean(ID_LCOMM, IDFLD, FXN, FYN, FAN)
Definition: w3ogcmmd.F90:407
w3adatmd::thm
real, dimension(:), pointer thm
Definition: w3adatmd.F90:587
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3servmd::w3s2xy
subroutine w3s2xy(NSEA, MSEA, MX, MY, S, MAPSF, XY)
Definition: w3servmd.F90:337
w3adatmd::t0m1
real, dimension(:), pointer t0m1
Definition: w3adatmd.F90:587
w3gdatmd::mapsta
integer, dimension(:,:), pointer mapsta
Definition: w3gdatmd.F90:1163
w3adatmd::tauocx
real, dimension(:), pointer tauocx
Definition: w3adatmd.F90:607
w3adatmd::ubd
real, dimension(:), pointer ubd
Definition: w3adatmd.F90:614