Module used for coupling applications between atmospheric model and WW3 with OASIS3-MCT.
More...
Module used for coupling applications between atmospheric model and WW3 with OASIS3-MCT.
- Author
- J. Pianezze
- Date
- Mar-2021
- Copyright
- Copyright 2009-2022 National Weather Service (NWS), National Oceanic and Atmospheric Administration. All rights reserved. WAVEWATCH III is a trademark of the NWS. No unauthorized use without permission.
| subroutine, public w3agcmmd::rcv_fields_from_atmos |
( |
integer, intent(in) |
ID_LCOMM, |
|
|
character(len=3), intent(in) |
IDFLD, |
|
|
real, dimension(:,:), intent(inout) |
FXN, |
|
|
real, dimension(:,:), intent(inout) |
FYN, |
|
|
real, dimension(:,:), intent(inout) |
FAN |
|
) |
| |
Receive coupling fields from atmospheric model.
- Parameters
-
| [in] | ID_LCOMM | MPI communicator. |
| [in] | IDFLD | Name of the exchange fields. |
| [in,out] | FXN | First exchange field. |
| [in,out] | FYN | Second exchange field. |
| [in,out] | FAN | Third exchange field. |
- Author
- J. Pianezze
- Date
- Mar-2021
Definition at line 235 of file w3agcmmd.F90.
296 INTEGER,
INTENT(IN) :: ID_LCOMM
297 CHARACTER(LEN=3),
INTENT(IN) :: IDFLD
298 REAL,
INTENT(INOUT) :: FXN(:,:), FYN(:,:), FAN(:,:)
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
312 rla_oasis_rcv(:,:) = 0.0
315 IF (idfld ==
'WND')
THEN
319 IF (
rcv_fld(ib_do)%CL_FIELD_NAME ==
'WW3__U10')
THEN
323 snd_buff(1:
nsea) = 0.0
326 snd_buff(ib_j) = tmp(ib_i)
329 CALL mpi_allreduce(snd_buff(1:
nsea), &
345 IF (
rcv_fld(ib_do)%CL_FIELD_NAME ==
'WW3__V10')
THEN
349 snd_buff(1:
nsea) = 0.0
352 snd_buff(ib_j) = tmp(ib_i)
355 CALL mpi_allreduce(snd_buff(1:
nsea), &
370 IF (idfld ==
'TAU')
THEN
374 IF (
rcv_fld(ib_do)%CL_FIELD_NAME ==
'WW3_UTAU')
THEN
378 snd_buff(1:
nsea) = 0.0
381 snd_buff(ib_j) = tmp(ib_i)
384 CALL mpi_allreduce(snd_buff(1:
nsea), &
400 IF (
rcv_fld(ib_do)%CL_FIELD_NAME ==
'WW3_VTAU')
THEN
404 snd_buff(1:
nsea) = 0.0
407 snd_buff(ib_j) = tmp(ib_i)
410 CALL mpi_allreduce(snd_buff(1:
nsea), &
425 IF (idfld ==
'RHO')
THEN
429 IF (
rcv_fld(ib_do)%CL_FIELD_NAME ==
'WW3_RHOA')
THEN
433 snd_buff(1:
nsea) = 0.0
436 snd_buff(ib_j) = tmp(ib_i)
439 CALL mpi_allreduce(snd_buff(1:
nsea), &
References w3oacpmd::cpl_oasis_rcv(), w3odatmd::iaproc, w3oacpmd::id_oasis_time, w3oacpmd::il_nb_rcv, w3gdatmd::mapsf, w3odatmd::naproc, w3gdatmd::nsea, w3gdatmd::nseal, w3gdatmd::nx, w3gdatmd::ny, w3oacpmd::rcv_fld, and w3servmd::w3s2xy().
Referenced by w3fldsmd::w3fldg().
| subroutine, public w3agcmmd::snd_fields_to_atmos |
Send coupling fields to atmospheric model.
- Author
- J. Pianezze
- Date
- Apr-2016
Definition at line 89 of file w3agcmmd.F90.
136 REAL(kind=8), dimension(
nseal,1) :: rla_oasis_snd
139 REAL(kind=8), dimension(
nseal) :: tmp
140 INTEGER :: JSEA, ISEA
149 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3_WSSU')
THEN
153 IF(
cx(isea) /= undef) tmp(jsea)=
cx(isea)
155 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
161 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3_WSSV')
THEN
165 IF(
cy(isea) /= undef) tmp(jsea)=
cy(isea)
167 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
173 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3_ACHA')
THEN
176 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
182 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3__AHS')
THEN
185 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
191 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3___FP')
THEN
194 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
200 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3___TP')
THEN
203 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
209 IF (
snd_fld(ib_do)%CL_FIELD_NAME ==
'WW3__FWS')
THEN
212 rla_oasis_snd(:,1) = dble(tmp(1:
nseal))
References w3adatmd::charn, w3oacpmd::cpl_oasis_snd(), w3adatmd::cx, w3adatmd::cy, w3adatmd::fp0, w3adatmd::hs, w3odatmd::iaproc, w3oacpmd::id_oasis_time, w3oacpmd::il_nb_snd, w3odatmd::naproc, w3gdatmd::nsea, w3gdatmd::nseal, w3oacpmd::snd_fld, and w3adatmd::tws.
Referenced by w3shel(), and w3wavemd::w3wave().