WAVEWATCH III  beta 0.0.1
w3igcmmd.F90
Go to the documentation of this file.
1 
6 
7 #include "w3macros.h"
8 
14 !/ ------------------------------------------------------------------- /
15 MODULE w3igcmmd
16  !/
17  !/ +-----------------------------------+
18  !/ | WAVEWATCH III NOAA/NCEP |
19  !/ | G. Boutin |
20  !/ | FORTRAN 90 |
21  !/ | Last update : Aug-2016 |
22  !/ +-----------------------------------+
23  !/
24  !/ Aug-2016 : Origination (G. Boutin) ( version 5.10 )
25  !/
26  !/ Copyright 2009-2012 National Weather Service (NWS),
27  !/ National Oceanic and Atmospheric Administration. All rights
28  !/ reserved. WAVEWATCH III is a trademark of the NWS.
29  !/ No unauthorized use without permission.
30  !/
31  ! 1. Purpose :
32  !
33  ! Module used for coupling applications between ice model and WW3 with OASIS3-MCT
34  !
35  ! 2. Variables and types :
36  !
37  ! 3. Subroutines and functions :
38  !
39  ! Name Type Scope Description
40  ! ----------------------------------------------------------------
41  ! SND_FIELDS_TO_ICE Subr. Public Send fields to ice model
42  ! RCV_FIELDS_FROM_ICE Subr. Public Receive fields from ice model
43  ! ----------------------------------------------------------------
44  !
45  ! 4. Subroutines and functions used :
46  !
47  ! Name Type Module Description
48  ! ----------------------------------------------------------------
49  ! CPL_OASIS_SEND Subr. W3OACPMD Send fields
50  ! CPL_OASIS_RECV Subr. W3OACPMD Receive fields
51  ! ----------------------------------------------------------------
52  !
53  ! 5. Remarks
54  ! 6. Switches :
55  ! 7. Source code :
56  !
57  !/ ------------------------------------------------------------------- /
58  !
59  IMPLICIT NONE
60  !
61  include "mpif.h"
62  !
63  PRIVATE
64  !
65  ! * Accessibility
66  PUBLIC snd_fields_to_ice
67  PUBLIC rcv_fields_from_ice
68  !
69 CONTAINS
70  !/ ------------------------------------------------------------------- /
76  SUBROUTINE snd_fields_to_ice()
77  !/
78  !/ +-----------------------------------+
79  !/ | WAVEWATCH III NOAA/NCEP |
80  !/ | G. Boutin |
81  !/ | FORTRAN 90 |
82  !/ | Last update : Aug-2016 |
83  !/ +-----------------------------------+
84  !/
85  !/ Aug-2016 : Origination (G. Boutin) ( version 5.10 )
86  !/
87  ! 1. Purpose :
88  !
89  ! Send coupling fields to ice model
90  !
91  ! 2. Method :
92  ! 3. Parameters :
93  ! 4. Subroutines used :
94  !
95  ! Name Type Module Description
96  ! -------------------------------------------------------------------
97  ! CPL_OASIS_SND Subr. W3OACPMD Send field to ice/ocean model
98  ! -------------------------------------------------------------------
99  !
100  ! 5. Called by :
101  !
102  ! Name Type Module Description
103  ! ------------------------------------------------------------------
104  ! W3WAVE Subr. W3WAVEMD Wave model
105  ! ------------------------------------------------------------------
106  !
107  ! 6. Error messages :
108  ! 7. Remarks :
109  ! 8. Structure :
110  ! 9. Switches :
111  ! 10. Source code :
112  !
113  !/ ------------------------------------------------------------------- /
114  !
116  USE w3gdatmd, ONLY: nseal, nsea
117  USE w3wdatmd, ONLY: icef
118  USE w3adatmd, ONLY: tauice
119  USE w3odatmd, ONLY: undef, naproc, iaproc
120  !
121  !/ ------------------------------------------------------------------- /
122  !/ Local parameters
123  !/
124  REAL(kind=8), dimension(nseal,1) :: rla_oasis_snd
125  INTEGER :: ib_do, ndso
126  LOGICAL :: ll_action
127  REAL(kind=8), dimension(nseal) :: tmp
128  INTEGER :: jsea, isea
129  !
130  !----------------------------------------------------------------------
131  ! * Executable part
132  !
133  DO ib_do = 1, il_nb_snd
134 
135  SELECT CASE(snd_fld(ib_do)%CL_FIELD_NAME)
136 
137  !
138  ! Ice floe diameters (m)
139  ! ---------------------------------------------------------------------
140  CASE ('WW3_ICEF')
141  tmp(1:nseal) = 0.0
142  DO jsea=1, nseal
143  isea=iaproc+(jsea-1)*naproc
144  IF(icef(isea) /= undef) tmp(jsea)=icef(isea)
145  END DO
146  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
147  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
148 
149  CASE ('WW3_TWIX')
150  tmp(1:nseal) = 0.0
151  WHERE(tauice(1:nseal,1) /= undef) tmp(1:nseal)=tauice(1:nseal,1)
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 
155  CASE ('WW3_TWIY')
156  tmp(1:nseal) = 0.0
157  WHERE(tauice(1:nseal,2) /= undef) tmp(1:nseal)=tauice(1:nseal,2)
158  rla_oasis_snd(:,1) = dble(tmp(1:nseal))
159  CALL cpl_oasis_snd(ib_do, id_oasis_time, rla_oasis_snd, ll_action)
160 
161  END SELECT
162 
163 
164  ENDDO
165  !
166  !/ ------------------------------------------------------------------- /
167  END SUBROUTINE snd_fields_to_ice
168  !/ ------------------------------------------------------------------- /
180  SUBROUTINE rcv_fields_from_ice(ID_LCOMM, IDFLD, FXN, FYN, FAN)
181  !/
182  !/ +-----------------------------------+
183  !/ | WAVEWATCH III NOAA/NCEP |
184  !/ | G. Boutin |
185  !/ | FORTRAN 90 |
186  !/ | Last update : April-2016 |
187  !/ +-----------------------------------+
188  !/
189  !/ Aug-2016 : Origination (G. Boutin) ( version 5.10 )
190  !/
191  ! 1. Purpose :
192  !
193  ! Receive coupling fields from ice model
194  !
195  ! 2. Method :
196  ! 3. Parameters :
197  !
198  ! Parameter list
199  ! ----------------------------------------------------------------
200  ! ID_LCOMM Char. I MPI communicator
201  ! IDFLD Int. I Name of the exchange fields
202  ! FXN Int. I/O First exchange field
203  ! FYN Int. I/O Second exchange field
204  ! FAN Int. I/O Third exchange field
205  ! ----------------------------------------------------------------
206  !
207  ! 4. Subroutines used :
208  !
209  ! Name Type Module Description
210  ! -------------------------------------------------------------------
211  ! CPL_OASIS_RCV Subr. W3OACPMD Receive fields from ice/ocean model
212  ! W3S2XY Subr. W3SERVMD Convert from storage (NSEA) to spatial grid (NX, NY)
213  ! -------------------------------------------------------------------
214  !
215  ! 5. Called by :
216  !
217  ! Name Type Module Description
218  ! ------------------------------------------------------------------
219  ! W3FLDG Subr. W3FLDSMD Manage input fields of depth,
220  ! current, wind and ice concentration
221  ! ------------------------------------------------------------------
222  !
223  ! 6. Error messages :
224  ! 7. Remarks :
225  ! 8. Structure :
226  ! 9. Switches :
227  ! 10. Source code :
228  !
229  !/ ------------------------------------------------------------------- /
230  !
232  USE w3gdatmd, ONLY: nx, ny, nseal, nsea, mapsf
233  USE w3odatmd, ONLY: naproc, iaproc
234  USE w3servmd, ONLY: w3s2xy
235  !
236  !/ ------------------------------------------------------------------- /
237  !/ Parameter list
238  !/
239  INTEGER, INTENT(IN) :: id_lcomm
240  CHARACTER(LEN=3), INTENT(IN) :: idfld
241  REAL, INTENT(INOUT) :: fxn(:,:), fyn(:,:), fan(:,:)
242  !
243  !/ ------------------------------------------------------------------- /
244  !/ Local parameters
245  !/
246  LOGICAL :: ll_action
247  INTEGER :: ib_do, ib_i, ib_j, il_err, ndso
248  REAL(kind=8), dimension(nseal,1) :: rla_oasis_rcv
249  REAL(kind=8), dimension(nseal) :: tmp
250  REAL, DIMENSION(1:NSEA) :: snd_buff,rcv_buff
251  !
252  !----------------------------------------------------------------------
253  ! * Executable part
254  !
255  rla_oasis_rcv(:,:) = 0.0
256  !
257  DO ib_do = 1, il_nb_rcv
258  IF (idfld == 'IC5') THEN
259  SELECT CASE (rcv_fld(ib_do)%CL_FIELD_NAME)
260  !
261  ! Ice floe diameters (m)
262  ! ----------------------------------------------------------------------
263  CASE ('WW3__IC5')
264 
265  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
266  IF (ll_action) THEN
267  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
268  snd_buff(1:nsea) = 0.0
269  DO ib_i = 1, nseal
270  ib_j = iaproc + (ib_i-1)*naproc
271  snd_buff(ib_j) = tmp(ib_i)
272  END DO
273  !
274  CALL mpi_allreduce(snd_buff(1:nsea), &
275  rcv_buff(1:nsea), &
276  nsea, &
277  mpi_real, &
278  mpi_sum, &
279  id_lcomm, &
280  il_err)
281  !
282  ! Convert from storage (NSEA) to spatial grid (NX, NY)
283  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fan)
284  !
285  END IF
286  END SELECT
287  !
288  ! Ice Concentration
289  ! ----------------------------------------------------------------------
290  ELSE IF (idfld == 'ICE') THEN
291  SELECT CASE (rcv_fld(ib_do)%CL_FIELD_NAME)
292  CASE ('WW3__ICE')
293  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
294  IF (ll_action) THEN
295  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
296  snd_buff(1:nsea) = 0.0
297  DO ib_i = 1, nseal
298  ib_j = iaproc + (ib_i-1)*naproc
299  snd_buff(ib_j) = tmp(ib_i)
300  END DO
301  !
302  !
303  CALL mpi_allreduce(snd_buff(1:nsea), &
304  rcv_buff(1:nsea), &
305  nsea, &
306  mpi_real, &
307  mpi_sum, &
308  id_lcomm, &
309  il_err)
310  !
311  ! Convert from storage (NSEA) to spatial grid (NX, NY)
312  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fan)
313  !
314  END IF
315  END SELECT
316  ! Ice Thickness
317  ! ----------------------------------------------------------------------
318  ELSE IF (idfld == 'IC1') THEN
319  SELECT CASE (rcv_fld(ib_do)%CL_FIELD_NAME)
320  CASE ('WW3__IC1')
321  CALL cpl_oasis_rcv(ib_do, id_oasis_time, rla_oasis_rcv, ll_action)
322  IF (ll_action) THEN
323  tmp(1:nseal) = rla_oasis_rcv(1:nseal,1)
324  snd_buff(1:nsea) = 0.0
325  DO ib_i = 1, nseal
326  ib_j = iaproc + (ib_i-1)*naproc
327  snd_buff(ib_j) = tmp(ib_i)
328  END DO
329 
330  CALL mpi_allreduce(snd_buff(1:nsea), &
331  rcv_buff(1:nsea), &
332  nsea, &
333  mpi_real, &
334  mpi_sum, &
335  id_lcomm, &
336  il_err)
337  !
338  ! Convert from storage (NSEA) to spatial grid (NX, NY)
339  CALL w3s2xy(nsea,nsea,nx,ny,rcv_buff(1:nsea),mapsf,fan)
340  ENDIF
341  !
342  END SELECT
343 
344  END IF
345  END DO
346  !
347 
348  !/ ------------------------------------------------------------------- /
349  END SUBROUTINE rcv_fields_from_ice
350  !/ ------------------------------------------------------------------- /
351  !/
352 END MODULE w3igcmmd
353 !/
354 !/ ------------------------------------------------------------------- /
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
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3wdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: w3wdatmd.F90:18
w3adatmd::tauice
real, dimension(:,:), pointer tauice
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
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
w3wdatmd::icef
real, dimension(:), pointer icef
Definition: w3wdatmd.F90:183
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3igcmmd::rcv_fields_from_ice
subroutine, public rcv_fields_from_ice(ID_LCOMM, IDFLD, FXN, FYN, FAN)
Receive coupling fields from ice model.
Definition: w3igcmmd.F90:181
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
w3oacpmd::il_nb_snd
integer, public il_nb_snd
Definition: w3oacpmd.F90:67
w3igcmmd::snd_fields_to_ice
subroutine, public snd_fields_to_ice()
Send coupling fields to ice model.
Definition: w3igcmmd.F90:77
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3odatmd
Definition: w3odatmd.F90:3
w3oacpmd
Definition: w3oacpmd.F90:3
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
w3igcmmd
Module used for coupling applications between ice model and WW3 with OASIS3-MCT.
Definition: w3igcmmd.F90:15
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
w3gdatmd
Definition: w3gdatmd.F90:16
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3servmd::w3s2xy
subroutine w3s2xy(NSEA, MSEA, MX, MY, S, MAPSF, XY)
Definition: w3servmd.F90:337