WAVEWATCH III  beta 0.0.1
w3agcmmd.F90
Go to the documentation of this file.
1 
8 
9 #include "w3macros.h"
10 !/ ------------------------------------------------------------------- /
23 MODULE w3agcmmd
24  !/
25  !/ +-----------------------------------+
26  !/ | WAVEWATCH III NOAA/NCEP |
27  !/ | J. Pianezze |
28  !/ | FORTRAN 90 |
29  !/ | Last update : Mar-2021 |
30  !/ +-----------------------------------+
31  !/
32  !/ Mar-2014 : Origination. ( version 4.18 )
33  !/ For upgrades see subroutines.
34  !/ Apr-2016 : Add comments (J. Pianezze) ( version 5.07 )
35  !/ Mar-2021 : Add TAUA and RHOA coupling ( version 7.13 )
36  !/
37  !/ Copyright 2009-2012 National Weather Service (NWS),
38  !/ National Oceanic and Atmospheric Administration. All rights
39  !/ reserved. WAVEWATCH III is a trademark of the NWS.
40  !/ No unauthorized use without permission.
41  !/
42  ! 1. Purpose :
43  !
44  ! Module used for coupling applications between atmospheric model and WW3 with OASIS3-MCT
45  !
46  ! 2. Variables and types :
47  !
48  ! 3. Subroutines and functions :
49  !
50  ! Name Type Scope Description
51  ! ----------------------------------------------------------------
52  ! SND_FIELDS_TO_ATMOS Subr. Public Send fields to atmos model
53  ! RCV_FIELDS_FROM_ATMOS Subr. Public Receive fields from atmos model
54  ! ----------------------------------------------------------------
55  !
56  ! 4. Subroutines and functions used :
57  !
58  ! Name Type Module Description
59  ! ----------------------------------------------------------------
60  ! CPL_OASIS_SEND Subr. W3OACPMD Send fields
61  ! CPL_OASIS_RECV Subr. W3OACPMD Receive fields
62  ! ----------------------------------------------------------------
63  !
64  ! 5. Remarks
65  ! 6. Switches :
66  ! 7. Source code :
67  !
68  !/ ------------------------------------------------------------------- /
69  !
70  IMPLICIT NONE
71  !
72  include "mpif.h"
73  !
74  PRIVATE
75  !
76  ! * Accessibility
77  PUBLIC snd_fields_to_atmos
79  !
80 CONTAINS
81  !/ ------------------------------------------------------------------- /
88  SUBROUTINE snd_fields_to_atmos()
89  !/
90  !/ +-----------------------------------+
91  !/ | WAVEWATCH III NOAA/NCEP |
92  !/ | J. Pianezze |
93  !/ | FORTRAN 90 |
94  !/ | Last update : Apr-2016 |
95  !/ +-----------------------------------+
96  !/
97  !/ Mar-2014 : Origination. ( version 4.18 )
98  !/ Apr-2016 : Add comments (J. Pianezze) ( version 5.07 )
99  !/
100  ! 1. Purpose :
101  !
102  ! Send coupling fields to atmospheric model
103  !
104  ! 2. Method :
105  ! 3. Parameters :
106  ! 4. Subroutines used :
107  !
108  ! Name Type Module Description
109  ! -------------------------------------------------------------------
110  ! CPL_OASIS_SND Subr. W3OACPMD Send field to atmos/ocean model
111  ! -------------------------------------------------------------------
112  !
113  ! 5. Called by :
114  !
115  ! Name Type Module Description
116  ! ------------------------------------------------------------------
117  ! W3WAVE Subr. W3WAVEMD Wave model
118  ! ------------------------------------------------------------------
119  !
120  ! 6. Error messages :
121  ! 7. Remarks :
122  ! 8. Structure :
123  ! 9. Switches :
124  ! 10. Source code :
125  !
126  !/ ------------------------------------------------------------------- /
127  !
129  USE w3gdatmd, ONLY: nseal, nsea
130  USE w3adatmd, ONLY: cx, cy, charn, hs, fp0, tws
131  USE w3odatmd, ONLY: undef, naproc, iaproc
132  !
133  !/ ------------------------------------------------------------------- /
134  !/ Local parameters
135  !/
136  REAL(kind=8), dimension(nseal,1) :: rla_oasis_snd
137  INTEGER :: ib_do
138  LOGICAL :: ll_action
139  REAL(kind=8), dimension(nseal) :: tmp
140  INTEGER :: jsea, isea
141  !
142  !----------------------------------------------------------------------
143  ! * Executable part
144  !
145  DO ib_do = 1, il_nb_snd
146  !
147  ! Ocean sea surface current (m.s-1) (u-component)
148  ! ---------------------------------------------------------------------
149  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_WSSU') THEN
150  tmp(1:nseal) = 0.0
151  DO jsea=1, nseal
152  isea=iaproc+(jsea-1)*naproc
153  IF(cx(isea) /= undef) tmp(jsea)=cx(isea)
154  END DO
155  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
156  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
157  ENDIF
158  !
159  ! Ocean sea surface current (m.s-1) (v-component)
160  ! ---------------------------------------------------------------------
161  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_WSSV') THEN
162  tmp(1:nseal) = 0.0
163  DO jsea=1, nseal
164  isea=iaproc+(jsea-1)*naproc
165  IF(cy(isea) /= undef) tmp(jsea)=cy(isea)
166  END DO
167  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
168  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
169  ENDIF
170  !
171  ! Charnock Coefficient (-)
172  ! ---------------------------------------------------------------------
173  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3_ACHA') THEN
174  tmp(1:nseal) = 0.0
175  WHERE(charn(1:nseal) /= undef) tmp(1:nseal)=charn(1:nseal)
176  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
177  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
178  ENDIF
179  !
180  ! Significant wave height (m)
181  ! ---------------------------------------------------------------------
182  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__AHS') THEN
183  tmp(1:nseal) = 0.0
184  WHERE(hs(1:nseal) /= undef) tmp(1:nseal)=hs(1:nseal)
185  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
186  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
187  ENDIF
188  !
189  ! Peak frequency (s-1)
190  ! ---------------------------------------------------------------------
191  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3___FP') THEN
192  tmp(1:nseal) = 0.0
193  WHERE(fp0(1:nseal) /= undef) tmp(1:nseal)=fp0(1:nseal)
194  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
195  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
196  ENDIF
197  !
198  ! Peak period (s)
199  ! ---------------------------------------------------------------------
200  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3___TP') THEN
201  tmp(1:nseal) = 0.0
202  WHERE(fp0(1:nseal) /= undef) tmp(1:nseal)=1./fp0(1:nseal)
203  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
204  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
205  ENDIF
206  !
207  ! Wind sea Mean period (s)
208  ! ---------------------------------------------------------------------
209  IF (snd_fld(ib_do)%CL_FIELD_NAME == 'WW3__FWS') THEN
210  tmp(1:nseal) = 0.0
211  WHERE(tws(1:nseal) /= undef) tmp(1:nseal)=tws(1:nseal)
212  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
213  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
214  ENDIF
215  !
216 
217  ENDDO
218  !
219  !/ ------------------------------------------------------------------- /
220  END SUBROUTINE snd_fields_to_atmos
233  !/ ------------------------------------------------------------------- /
234  SUBROUTINE rcv_fields_from_atmos(ID_LCOMM, IDFLD, FXN, FYN, FAN)
235  !/
236  !/ +-----------------------------------+
237  !/ | WAVEWATCH III NOAA/NCEP |
238  !/ | J. Pianezze |
239  !/ | FORTRAN 90 |
240  !/ | Last update : Mar-2021 |
241  !/ +-----------------------------------+
242  !/
243  !/ Mar-2014 : Origination. ( version 4.18 )
244  !/ Apr-2015 : Modification (M. Accensi) ( version 5.07 )
245  !/ Apr-2016 : Add comments (J. Pianezze) ( version 5.07 )
246  !/ Mar-2021 : Add TAUA and RHOA coupling ( version 7.13 )
247  !/
248  ! 1. Purpose :
249  !
250  ! Receive coupling fields from atmospheric model
251  !
252  ! 2. Method :
253  ! 3. Parameters :
254  !
255  ! Parameter list
256  ! ----------------------------------------------------------------
257  ! ID_LCOMM Char. I MPI communicator
258  ! IDFLD Int. I Name of the exchange fields
259  ! FXN Int. I/O First exchange field
260  ! FYN Int. I/O Second exchange field
261  ! FAN Int. I/O Third exchange field
262  ! ----------------------------------------------------------------
263  !
264  ! 4. Subroutines used :
265  !
266  ! Name Type Module Description
267  ! -------------------------------------------------------------------
268  ! CPL_OASIS_RCV Subr. W3OACPMD Receive fields from atmos/ocean model
269  ! W3S2XY Subr. W3SERVMD Convert from storage (NSEA) to spatial grid (NX, NY)
270  ! -------------------------------------------------------------------
271  !
272  ! 5. Called by :
273  !
274  ! Name Type Module Description
275  ! ------------------------------------------------------------------
276  ! W3FLDG Subr. W3FLDSMD Manage input fields of depth,
277  ! current, wind and ice concentration
278  ! ------------------------------------------------------------------
279  !
280  ! 6. Error messages :
281  ! 7. Remarks :
282  ! 8. Structure :
283  ! 9. Switches :
284  ! 10. Source code :
285  !
286  !/ ------------------------------------------------------------------- /
287  !
289  USE w3gdatmd, ONLY: nx, ny, nseal, nsea, mapsf
290  USE w3odatmd, ONLY: naproc, iaproc
291  USE w3servmd, ONLY: w3s2xy
292  !
293  !/ ------------------------------------------------------------------- /
294  !/ Parameter list
295  !/
296  INTEGER, INTENT(IN) :: id_lcomm
297  CHARACTER(LEN=3), INTENT(IN) :: idfld
298  REAL, INTENT(INOUT) :: fxn(:,:), fyn(:,:), fan(:,:)
299  !
300  !/ ------------------------------------------------------------------- /
301  !/ Local parameters
302  !/
303  LOGICAL :: ll_action
304  INTEGER :: ib_do, ib_i, ib_j, il_err
305  REAL(kind=8), dimension(nseal,1) :: rla_oasis_rcv
306  REAL(kind=8), dimension(nseal) :: tmp
307  REAL, DIMENSION(1:NSEA) :: snd_buff,rcv_buff
308  !
309  !----------------------------------------------------------------------
310  ! * Executable part
311  !
312  rla_oasis_rcv(:,:) = 0.0
313  !
314  DO ib_do = 1, il_nb_rcv
315  IF (idfld == 'WND') THEN
316  !
317  ! Wind speed at 10m (m.s-1) (u-component)
318  ! ----------------------------------------------------------------------
319  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3__U10') THEN
320  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
321  IF (ll_action) THEN
322  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
323  snd_buff(1:nsea) = 0.0
324  DO ib_i = 1, nseal
325  ib_j = iaproc + (ib_i-1)*naproc
326  snd_buff(ib_j) = tmp(ib_i)
327  ENDDO
328  !
329  CALL mpi_allreduce(snd_buff(1:nsea), &
330  rcv_buff(1:nsea), &
331  nsea, &
332  mpi_real, &
333  mpi_sum, &
334  id_lcomm, &
335  il_err)
336  !
337  ! Convert from storage (NSEA) to spatial grid (NX, NY)
338  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fxn)
339  !
340  ENDIF
341  ENDIF
342  !
343  ! Wind speed at 10m (m.s-1) (v-component)
344  ! ----------------------------------------------------------------------
345  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3__V10') THEN
346  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
347  IF (ll_action) THEN
348  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
349  snd_buff(1:nsea) = 0.0
350  DO ib_i = 1, nseal
351  ib_j = iaproc + (ib_i-1)*naproc
352  snd_buff(ib_j) = tmp(ib_i)
353  END DO
354  !
355  CALL mpi_allreduce(snd_buff(1:nsea), &
356  rcv_buff(1:nsea), &
357  nsea, &
358  mpi_real, &
359  mpi_sum, &
360  id_lcomm, &
361  il_err)
362  !
363  ! Convert from storage (NSEA) to spatial grid (NX, NY)
364  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fyn)
365  !
366  ENDIF
367  ENDIF
368  !
369  ENDIF
370  IF (idfld == 'TAU') THEN
371  !
372  ! Atmospheric momentum (Pa) (u-component)
373  ! ----------------------------------------------------------------------
374  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_UTAU') THEN
375  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
376  IF (ll_action) THEN
377  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
378  snd_buff(1:nsea) = 0.0
379  DO ib_i = 1, nseal
380  ib_j = iaproc + (ib_i-1)*naproc
381  snd_buff(ib_j) = tmp(ib_i)
382  ENDDO
383  !
384  CALL mpi_allreduce(snd_buff(1:nsea), &
385  rcv_buff(1:nsea), &
386  nsea, &
387  mpi_real, &
388  mpi_sum, &
389  id_lcomm, &
390  il_err)
391  !
392  ! Convert from storage (NSEA) to spatial grid (NX, NY)
393  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fxn)
394  !
395  ENDIF
396  ENDIF
397  !
398  ! Atmospheric momentum (Pa) (v-component)
399  ! ----------------------------------------------------------------------
400  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_VTAU') THEN
401  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
402  IF (ll_action) THEN
403  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
404  snd_buff(1:nsea) = 0.0
405  DO ib_i = 1, nseal
406  ib_j = iaproc + (ib_i-1)*naproc
407  snd_buff(ib_j) = tmp(ib_i)
408  END DO
409  !
410  CALL mpi_allreduce(snd_buff(1:nsea), &
411  rcv_buff(1:nsea), &
412  nsea, &
413  mpi_real, &
414  mpi_sum, &
415  id_lcomm, &
416  il_err)
417  !
418  ! Convert from storage (NSEA) to spatial grid (NX, NY)
419  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fyn)
420  !
421  ENDIF
422  ENDIF
423  !
424  ENDIF
425  IF (idfld == 'RHO') THEN
426  !
427  ! Air density (kg.m-3)
428  ! ----------------------------------------------------------------------
429  IF (rcv_fld(ib_do)%CL_FIELD_NAME == 'WW3_RHOA') THEN
430  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
431  IF (ll_action) THEN
432  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
433  snd_buff(1:nsea) = 0.0
434  DO ib_i = 1, nseal
435  ib_j = iaproc + (ib_i-1)*naproc
436  snd_buff(ib_j) = tmp(ib_i)
437  ENDDO
438  !
439  CALL mpi_allreduce(snd_buff(1:nsea), &
440  rcv_buff(1:nsea), &
441  nsea, &
442  mpi_real, &
443  mpi_sum, &
444  id_lcomm, &
445  il_err)
446  !
447  ! Convert from storage (NSEA) to spatial grid (NX, NY)
448  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fan)
449  !
450  ENDIF
451  ENDIF
452  ENDIF
453  ENDDO
454  !/ ------------------------------------------------------------------- /
455  END SUBROUTINE rcv_fields_from_atmos
456  !/ ------------------------------------------------------------------- /
457  !/
458 END MODULE w3agcmmd
459 !/
460 !/ ------------------------------------------------------------------- /
w3gdatmd::nseal
integer, pointer nseal
Definition: w3gdatmd.F90:1097
include
cmake src_list cmake include(${CMAKE_CURRENT_SOURCE_DIR}/cmake/check_switches.cmake) check_switches("$
Definition: CMakeLists.txt:15
w3agcmmd
Module used for coupling applications between atmospheric model and WW3 with OASIS3-MCT.
Definition: w3agcmmd.F90:23
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::tws
real, dimension(:), pointer tws
Definition: w3adatmd.F90:603
w3oacpmd::cpl_oasis_snd
subroutine, public cpl_oasis_snd(ID_NB, ID_TIME, RDA_FIELD, LD_ACTION)
Definition: w3oacpmd.F90:563
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
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3agcmmd::rcv_fields_from_atmos
subroutine, public rcv_fields_from_atmos(ID_LCOMM, IDFLD, FXN, FYN, FAN)
Receive coupling fields from atmospheric model.
Definition: w3agcmmd.F90:235
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
w3agcmmd::snd_fields_to_atmos
subroutine, public snd_fields_to_atmos()
Send coupling fields to atmospheric model.
Definition: w3agcmmd.F90:89
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3odatmd
Definition: w3odatmd.F90:3
w3oacpmd
Definition: w3oacpmd.F90:3
w3adatmd::cy
real, dimension(:), pointer cy
Definition: w3adatmd.F90:584
w3gdatmd::mapsf
integer, dimension(:,:), pointer mapsf
Definition: w3gdatmd.F90:1163
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
w3oacpmd::rcv_fld
type(cpl_field), dimension(ip_maxfld), public rcv_fld
Definition: w3oacpmd.F90:76
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
w3adatmd::fp0
real, dimension(:), pointer fp0
Definition: w3adatmd.F90:587
w3gdatmd
Definition: w3gdatmd.F90:16
w3adatmd::cx
real, dimension(:), pointer cx
Definition: w3adatmd.F90:584
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3servmd::w3s2xy
subroutine w3s2xy(NSEA, MSEA, MX, MY, S, MAPSF, XY)
Definition: w3servmd.F90:337