WAVEWATCH III  beta 0.0.1
w3wavset Module Reference

Implicit solution of wave setup problem following Dingemans for structured and unstructured grids. More...

Functions/Subroutines

subroutine differentiate_xydir_native (VAR, DVDX, DVDY)
 Differentiate xy, using linear shape function. More...
 
subroutine differentiate_xydir_mapsta (VAR, DVDX, DVDY)
 Differentiate xy based on mapsta, using linear shape function. More...
 
subroutine differentiate_xydir (VAR, DVDX, DVDY)
 Driver routine for xydir. More...
 
subroutine trig_compute_lh_stress (F_X, F_Y, DWNX)
 Setup boundary pointer. More...
 
subroutine trig_compute_diff (IE, I1, UGRAD, VGRAD)
 Differentiate other way around. More...
 
subroutine trig_wave_setup_compute_system (ASPAR, B, FX, FY, DWNX, ACTIVE, ACTIVESEC)
 Setup system matrix for solutions of wave setup eq. More...
 
subroutine trig_wave_setup_apply_precond (ASPAR, TheIn, TheOut, ACTIVE, ACTIVESEC)
 Preconditioner. More...
 
subroutine trig_wave_setup_apply_fct (ASPAR, TheIn, TheOut, ACTIVE, ACTIVESEC)
 
subroutine trig_wave_setup_scalar_prod (V1, V2, eScal)
 Scalar product plus exchange. More...
 
subroutine trig_wave_setup_solve_poisson_neumann_dir (ASPAR, B, TheOut, ACTIVE, ACTIVESEC)
 Poisson equation solver. More...
 
subroutine trig_set_meanvalue_to_zero (TheVar)
 Set mean value. More...
 
subroutine compute_active_node (DWNX, ACTIVE)
 Compute active node for setup comp. More...
 
subroutine trig_wave_setup_computation
 Setup computation. More...
 
subroutine preparation_fd_scheme (IMOD)
 Wave setup for FD grids. More...
 
subroutine fd_wave_setup_apply_fct (ASPAR, TheIn, TheOut)
 Compute off diagonal for FD grids. More...
 
subroutine fd_wave_setup_apply_precond (ASPAR, TheIn, TheOut)
 Preconditioning for FD grids. More...
 
subroutine fd_collect_sxx_xy_yy (SXX_t, SXY_t, SYY_t)
 Radiation stresses for FD grids. More...
 
subroutine fd_compute_lh_stress (SXX_t, SXY_t, SYY_t, FX, FY)
 Setup fluxes. More...
 
subroutine fd_compute_diff (IEDGE, ISEA, UGRAD, VGRAD, dist)
 Differences on FD grids. More...
 
subroutine fd_wave_setup_compute_system (ASPAR, B, FX, FY)
 Setup matrix on FD grids. More...
 
subroutine fd_wave_setup_scalar_prod (V1, V2, eScal)
 Scalar product. More...
 
subroutine fd_wave_setup_solve_poisson_neumann_dir (ASPAR, B, TheOut)
 Poisson solver on FD grids. More...
 
subroutine fd_set_meanvalue_to_zero (TheVar)
 Set mean value. More...
 
subroutine fd_wave_setup_computation
 Wave setup comp on FD grids. More...
 
subroutine wave_setup_computation
 General driver. More...
 

Variables

integer, save ient = 0
 
logical do_wave_setup = .TRUE.
 

Detailed Description

Implicit solution of wave setup problem following Dingemans for structured and unstructured grids.

Author
Aron Roland
Mathieu Dutour-Sikiric
Date
1-Jun-2018

Function/Subroutine Documentation

◆ compute_active_node()

subroutine w3wavset::compute_active_node ( real(rkind), dimension(npa), intent(in)  DWNX,
integer, dimension(npa), intent(out)  ACTIVE 
)

Compute active node for setup comp.

Parameters
[in]DWNX
[out]ACTIVE
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
1-May-2018

Definition at line 1494 of file w3wavset.F90.

1494 !/
1495 !/ +-----------------------------------+
1496 !/ | WAVEWATCH III NOAA/NCEP |
1497 !/ | |
1498 !/ | Aron Roland (BGS IT&E GmbH) |
1499 !/ | Mathieu Dutour-Sikiric (IRB) |
1500 !/ | |
1501 !/ | FORTRAN 90 |
1502 !/ | Last update : 01-Mai-2018 |
1503 !/ +-----------------------------------+
1504 !/
1505 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1506 !/
1507 ! 1. Purpose : Compute active node for setup comp.
1508 ! 2. Method :
1509 ! 3. Parameters :
1510 !
1511 ! Parameter list
1512 ! ----------------------------------------------------------------
1513 ! ----------------------------------------------------------------
1514 !
1515 ! 4. Subroutines used :
1516 !
1517 ! Name Type Module Description
1518 ! ----------------------------------------------------------------
1519 ! STRACE Subr. W3SERVMD Subroutine tracing.
1520 ! ----------------------------------------------------------------
1521 !
1522 ! 5. Called by :
1523 !
1524 ! Name Type Module Description
1525 ! ----------------------------------------------------------------
1526 ! ----------------------------------------------------------------
1527 !
1528 ! 6. Error messages :
1529 ! 7. Remarks
1530 ! 8. Structure :
1531 ! 9. Switches :
1532 !
1533 ! !/S Enable subroutine tracing.
1534 !
1535 ! 10. Source code :
1536 !
1537 !/ ------------------------------------------------------------------- /
1538 #ifdef W3_S
1539  USE w3servmd, ONLY: strace
1540 #endif
1541 !
1542  USE w3gdatmd, ONLY : crit_dep_stp
1543  USE yownodepool, only: pdlib_nnz, pdlib_ia, pdlib_ja, iplg, npa, np
1544  USE w3odatmd, only : iaproc
1545  IMPLICIT NONE
1546 !/
1547 !/ ------------------------------------------------------------------- /
1548 !/ Parameter list
1549 !/
1550 !/ ------------------------------------------------------------------- /
1551 !/ Local PARAMETERs
1552 !/
1553 #ifdef W3_S
1554  INTEGER, SAVE :: IENT = 0
1555 #endif
1556 !/
1557 !/ ------------------------------------------------------------------- /
1558 !/
1559  REAL(rkind), INTENT(in) :: DWNX(npa)
1560  INTEGER, INTENT(out) :: ACTIVE(npa)
1561  INTEGER IP, eAct
1562 #ifdef W3_DEBUGSTP
1563  INTEGER nbActive
1564 #endif
1565 #ifdef W3_S
1566  CALL strace (ient, 'VA_SETUP_IOBPD')
1567 #endif
1568 #ifdef W3_DEBUGSTP
1569  nbactive=0
1570 #endif
1571  DO ip=1,npa
1572  IF (dwnx(ip) .ge. crit_dep_stp) THEN
1573  eact=1
1574  ELSE
1575  eact=0
1576  END IF
1577 #ifdef W3_DEBUGSTP
1578  nbactive=nbactive + eact
1579 #endif
1580  active(ip)=eact
1581  END DO
1582 #ifdef W3_DEBUGSTP
1583  WRITE(740+iaproc,*) 'min/max(DWNX)=', minval(dwnx), maxval(dwnx)
1584  WRITE(740+iaproc,*) 'CRIT_DEP_STP=', crit_dep_stp
1585  WRITE(740+iaproc,*) 'nbActive=', nbactive, ' npa=', npa
1586  FLUSH(740+iaproc)
1587 #endif

References w3gdatmd::crit_dep_stp, w3odatmd::iaproc, yownodepool::iplg, yownodepool::np, yownodepool::npa, yownodepool::pdlib_ia, yownodepool::pdlib_ja, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by trig_wave_setup_computation().

◆ differentiate_xydir()

subroutine w3wavset::differentiate_xydir ( real(rkind), dimension(npa), intent(in)  VAR,
real(rkind), dimension(npa), intent(out)  DVDX,
real(rkind), dimension(npa), intent(out)  DVDY 
)

Driver routine for xydir.

Parameters
[in]VAR
[out]DVDX
[out]DVDY
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 364 of file w3wavset.F90.

364 !/
365 !/ +-----------------------------------+
366 !/ | WAVEWATCH III NOAA/NCEP |
367 !/ | |
368 !/ | Mathieu Dutour-Sikiric (IRB) |
369 !/ | Aron Roland (BGS IT&E GmbH) |
370 !/ | FORTRAN 90 |
371 !/ | Last update : 01-Mai-2018 |
372 !/ +-----------------------------------+
373 !/
374 !/ 01-Mai-2018 : Origination. ( version 6.04 )
375 !/
376 ! 1. Purpose : Driver routine for xydir
377 ! 2. Method :
378 ! 3. Parameters :
379 !
380 ! Parameter list
381 ! ----------------------------------------------------------------
382 ! ----------------------------------------------------------------
383 !
384 ! 4. Subroutines used :
385 !
386 ! Name Type Module Description
387 ! ----------------------------------------------------------------
388 ! STRACE Subr. W3SERVMD Subroutine tracing.
389 ! ----------------------------------------------------------------
390 !
391 ! 5. Called by :
392 !
393 ! Name Type Module Description
394 ! ----------------------------------------------------------------
395 ! ----------------------------------------------------------------
396 !
397 ! 6. Error messages :
398 ! 7. Remarks
399 ! 8. Structure :
400 ! 9. Switches :
401 !
402 ! !/S Enable subroutine tracing.
403 !
404 ! 10. Source code :
405 !
406 !/ ------------------------------------------------------------------- /
407 #ifdef W3_S
408  USE w3servmd, ONLY: strace
409 #endif
410 !
411  use yownodepool, only: npa
412  IMPLICIT NONE
413 !/
414 !/ ------------------------------------------------------------------- /
415 !/ Parameter list
416 !/
417 !/ ------------------------------------------------------------------- /
418 !/ Local PARAMETERs
419 !/
420 #ifdef W3_S
421  INTEGER, SAVE :: IENT = 0
422 #endif
423  REAL(rkind), INTENT(IN) :: VAR(npa)
424  REAL(rkind), INTENT(OUT) :: DVDX(npa), DVDY(npa)
425 !/
426 !/ ------------------------------------------------------------------- /
427 !/
428 #ifdef W3_S
429  CALL strace (ient, 'VA_SETUP_IOBPD')
430 #endif
431 !
432 
433  CALL differentiate_xydir_mapsta(var, dvdx, dvdy)
434 ! CALL DIFFERENTIATE_XYDIR_NATIVE(VAR, DVDX, DVDY)

References differentiate_xydir_mapsta(), yownodepool::npa, and w3servmd::strace().

Referenced by trig_compute_lh_stress().

◆ differentiate_xydir_mapsta()

subroutine w3wavset::differentiate_xydir_mapsta ( real(rkind), dimension(npa), intent(in)  VAR,
real(rkind), dimension(npa), intent(out)  DVDX,
real(rkind), dimension(npa), intent(out)  DVDY 
)

Differentiate xy based on mapsta, using linear shape function.

Parameters
[in]VAR
[out]DVDX
[out]DVDY
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
1-May-2018

Definition at line 226 of file w3wavset.F90.

226 !/
227 !/ +-----------------------------------+
228 !/ | WAVEWATCH III NOAA/NCEP |
229 !/ | |
230 !/ | Aron Roland (BGS IT&E GmbH) |
231 !/ | Mathieu Dutour-Sikiric (IRB) |
232 !/ | |
233 !/ | FORTRAN 90 |
234 !/ | Last update : 01-Mai-2018 |
235 !/ +-----------------------------------+
236 !/
237 !/ 01-Mai-2018 : Origination. ( version 6.04 )
238 !/
239 ! 1. Purpose : differentiate xy based on mapsta
240 ! 2. Method : linear shape function
241 ! 3. Parameters :
242 !
243 ! Parameter list
244 ! ----------------------------------------------------------------
245 ! ----------------------------------------------------------------
246 !
247 ! 4. Subroutines used :
248 !
249 ! Name Type Module Description
250 ! ----------------------------------------------------------------
251 ! STRACE Subr. W3SERVMD Subroutine tracing.
252 ! ----------------------------------------------------------------
253 !
254 ! 5. Called by :
255 !
256 ! Name Type Module Description
257 ! ----------------------------------------------------------------
258 ! ----------------------------------------------------------------
259 !
260 ! 6. Error messages :
261 ! 7. Remarks
262 ! 8. Structure :
263 ! 9. Switches :
264 !
265 ! !/S Enable subroutine tracing.
266 !
267 ! 10. Source code :
268 !
269 !/ ------------------------------------------------------------------- /
270 #ifdef W3_S
271  USE w3servmd, ONLY: strace
272 #endif
273 !
275  use yownodepool, only : pdlib_ien, pdlib_tria, npa, iplg
276  use yowelementpool, only : ine, ne
277  USE w3gdatmd, ONLY : mapsta
278  USE w3parall, only: init_get_isea
279  IMPLICIT NONE
280 !/
281 !/ ------------------------------------------------------------------- /
282 !/ Parameter list
283 !/
284 !/ ------------------------------------------------------------------- /
285 !/ Local PARAMETERs
286 !/
287 #ifdef W3_S
288  INTEGER, SAVE :: IENT = 0
289 #endif
290 !/
291 !/ ------------------------------------------------------------------- /
292 !/
293  REAL(rkind), INTENT(IN) :: VAR(npa)
294  REAL(rkind), INTENT(OUT) :: DVDX(npa), DVDY(npa)
295  INTEGER :: NI(3)
296  INTEGER :: IE, I1, I2, I3, IP, IX
297  REAL(rkind) :: DEDY(3),DEDX(3)
298  REAL(rkind) :: DVDXIE, DVDYIE
299  REAL(rkind) :: WEI(npa), eW
300  INTEGER :: IX1, IX2, IX3, ISEA
301 #ifdef W3_S
302  CALL strace (ient, 'VA_SETUP_IOBPD')
303 #endif
304  wei = 0.0
305  dvdx = 0.0
306  dvdy = 0.0
307 
308  DO ie = 1, ne
309  ni = ine(:,ie)
310  i1 = ine(1,ie)
311  i2 = ine(2,ie)
312  i3 = ine(3,ie)
313  ix1=iplg(i1)
314  ix2=iplg(i2)
315  ix3=iplg(i3)
316  IF ((mapsta(1,ix1) .gt. 0).and.(mapsta(1,ix2) .gt. 0).and.(mapsta(1,ix3) .gt. 0)) THEN
317  wei(ni) = wei(ni) + 2.*pdlib_tria(ie)
318  dedx(1) = pdlib_ien(1,ie)
319  dedx(2) = pdlib_ien(3,ie)
320  dedx(3) = pdlib_ien(5,ie)
321  dedy(1) = pdlib_ien(2,ie)
322  dedy(2) = pdlib_ien(4,ie)
323  dedy(3) = pdlib_ien(6,ie)
324  dvdxie = dot_product( var(ni),dedx)
325  dvdyie = dot_product( var(ni),dedy)
326  dvdx(ni) = dvdx(ni) + dvdxie
327  dvdy(ni) = dvdy(ni) + dvdyie
328  END IF
329  END DO
330  DO ip=1,npa
331  ix=iplg(ip)
332  ew=wei(ip)
333  IF (ew .gt. 0 .and. mapsta(1,ix) .gt. 0) THEN
334  dvdx(ip)=dvdx(ip) / ew
335  dvdy(ip)=dvdy(ip) / ew
336  ELSE
337  dvdx(ip)=0.
338  dvdy(ip)=0.
339  ENDIF
340  END DO
341  DO ip=1,npa
342  ix=iplg(ip)
343  IF (mapsta(1,ix) .lt. 0) THEN
344  dvdx(ip)=0.
345  dvdy(ip)=0.
346  END IF
347  END DO
348  CALL pdlib_exchange1dreal(dvdx)
349  CALL pdlib_exchange1dreal(dvdy)

References yowelementpool::ine, w3parall::init_get_isea(), yownodepool::iplg, w3gdatmd::mapsta, yowelementpool::ne, yownodepool::npa, yowexchangemodule::pdlib_exchange1dreal(), yownodepool::pdlib_ien, yownodepool::pdlib_tria, and w3servmd::strace().

Referenced by differentiate_xydir().

◆ differentiate_xydir_native()

subroutine w3wavset::differentiate_xydir_native ( real(rkind), dimension(npa), intent(in)  VAR,
real(rkind), dimension(npa), intent(out)  DVDX,
real(rkind), dimension(npa), intent(out)  DVDY 
)

Differentiate xy, using linear shape function.

Parameters
[in]VAR
[out]DVDX
[out]DVDY
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
1-May-2018

Definition at line 106 of file w3wavset.F90.

106 !/
107 !/ +-----------------------------------+
108 !/ | WAVEWATCH III NOAA/NCEP |
109 !/ | |
110 !/ | Aron Roland (BGS IT&E GmbH) |
111 !/ | Mathieu Dutour-Sikiric (IRB) |
112 !/ | |
113 !/ | FORTRAN 90 |
114 !/ | Last update : 01-Mai-2018 |
115 !/ +-----------------------------------+
116 !/
117 !/ 01-Mai-2018 : Origination. ( version 6.04 )
118 !/
119 ! 1. Purpose : differentiate xy
120 ! 2. Method : linear shape function
121 ! 3. Parameters :
122 !
123 ! Parameter list
124 ! ----------------------------------------------------------------
125 ! ----------------------------------------------------------------
126 !
127 ! 4. Subroutines used :
128 !
129 ! Name Type Module Description
130 ! ----------------------------------------------------------------
131 ! STRACE Subr. W3SERVMD Subroutine tracing.
132 ! ----------------------------------------------------------------
133 !
134 ! 5. Called by :
135 !
136 ! Name Type Module Description
137 ! ----------------------------------------------------------------
138 ! ----------------------------------------------------------------
139 !
140 ! 6. Error messages :
141 ! 7. Remarks
142 ! 8. Structure :
143 ! 9. Switches :
144 !
145 ! !/S Enable subroutine tracing.
146 !
147 ! 10. Source code :
148 !
149 !/ ------------------------------------------------------------------- /
150 #ifdef W3_S
151  USE w3servmd, ONLY: strace
152 #endif
153 !
155  use yownodepool, only : pdlib_ien, pdlib_tria, npa
156  use yowelementpool, only : ine, ne
157  USE w3gdatmd, ONLY : mapsta
158  IMPLICIT NONE
159 !/
160 !/ ------------------------------------------------------------------- /
161 !/ Parameter list
162 !/
163 !/ ------------------------------------------------------------------- /
164 !/ Local PARAMETERs
165 !/
166 #ifdef W3_S
167  INTEGER, SAVE :: IENT = 0
168 #endif
169 !/
170 !/ ------------------------------------------------------------------- /
171 !/
172 !
173  REAL(rkind), INTENT(IN) :: VAR(npa)
174  REAL(rkind), INTENT(OUT) :: DVDX(npa), DVDY(npa)
175  INTEGER :: NI(3)
176  INTEGER :: IE, I1, I2, I3, IP
177  REAL(rkind) :: DEDY(3),DEDX(3)
178  REAL(rkind) :: DVDXIE, DVDYIE
179  REAL(rkind) :: WEI(npa), eW
180  INTEGER :: IX
181 #ifdef W3_S
182  CALL strace (ient, 'VA_SETUP_IOBPD')
183 #endif
184  wei = 0.0
185  dvdx = 0.0
186  dvdy = 0.0
187 
188  DO ie = 1, ne
189  ni = ine(:,ie)
190  i1 = ine(1,ie)
191  i2 = ine(2,ie)
192  i3 = ine(3,ie)
193  wei(ni) = wei(ni) + 2.*pdlib_tria(ie)
194  dedx(1) = pdlib_ien(1,ie)
195  dedx(2) = pdlib_ien(3,ie)
196  dedx(3) = pdlib_ien(5,ie)
197  dedy(1) = pdlib_ien(2,ie)
198  dedy(2) = pdlib_ien(4,ie)
199  dedy(3) = pdlib_ien(6,ie)
200  dvdxie = dot_product( var(ni),dedx)
201  dvdyie = dot_product( var(ni),dedy)
202  dvdx(ni) = dvdx(ni) + dvdxie
203  dvdy(ni) = dvdy(ni) + dvdyie
204  END DO
205  DO ix=1,npa
206  ew=wei(ix)
207  dvdx(ix)=dvdx(ix) / ew
208  dvdy(ix)=dvdy(ix) / ew
209  END DO
210  CALL pdlib_exchange1dreal(dvdx)
211  CALL pdlib_exchange1dreal(dvdy)

References yowelementpool::ine, w3gdatmd::mapsta, yowelementpool::ne, yownodepool::npa, yowexchangemodule::pdlib_exchange1dreal(), yownodepool::pdlib_ien, yownodepool::pdlib_tria, and w3servmd::strace().

◆ fd_collect_sxx_xy_yy()

subroutine w3wavset::fd_collect_sxx_xy_yy ( real(rkind), dimension(nsea), intent(out)  SXX_t,
real(rkind), dimension(nsea), intent(out)  SXY_t,
real(rkind), dimension(nsea), intent(out)  SYY_t 
)

Radiation stresses for FD grids.

Parameters
[out]SXX_t
[out]SXY_t
[out]SYY_t
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2153 of file w3wavset.F90.

2153 !/
2154 !/ +-----------------------------------+
2155 !/ | WAVEWATCH III NOAA/NCEP |
2156 !/ | |
2157 !/ | Mathieu Dutour-Sikiric (IRB) |
2158 !/ | Aron Roland (BGS IT&E GmbH) |
2159 !/ | |
2160 !/ | FORTRAN 90 |
2161 !/ | Last update : 01-Mai-2018 |
2162 !/ +-----------------------------------+
2163 !/
2164 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2165 !/
2166 ! 1. Purpose : Rad. stresses for FD grids
2167 ! 2. Method :
2168 ! 3. Parameters :
2169 !
2170 ! Parameter list
2171 ! ----------------------------------------------------------------
2172 ! ----------------------------------------------------------------
2173 !
2174 ! 4. Subroutines used :
2175 !
2176 ! Name Type Module Description
2177 ! ----------------------------------------------------------------
2178 ! STRACE Subr. W3SERVMD Subroutine tracing.
2179 ! ----------------------------------------------------------------
2180 !
2181 ! 5. Called by :
2182 !
2183 ! Name Type Module Description
2184 ! ----------------------------------------------------------------
2185 ! ----------------------------------------------------------------
2186 !
2187 ! 6. Error messages :
2188 ! 7. Remarks
2189 ! 8. Structure :
2190 ! 9. Switches :
2191 !
2192 ! !/S Enable subroutine tracing.
2193 !
2194 ! 10. Source code :
2195 !
2196 !/ ------------------------------------------------------------------- /
2197 #ifdef W3_S
2198  USE w3servmd, ONLY: strace
2199 #endif
2200 !
2201  USE w3adatmd, ONLY: sxx, sxy, syy
2202  USE w3gdatmd, ONLY: nsea, nseal
2203  USE w3odatmd, only : iaproc, naproc
2204  use yowdatapool, only: rtype, istatus
2205  USE w3adatmd, ONLY: mpi_comm_wcmp
2206  IMPLICIT NONE
2207 !/
2208 !/ ------------------------------------------------------------------- /
2209 !/ Parameter list
2210 !/
2211 !/ ------------------------------------------------------------------- /
2212 !/ Local PARAMETERs
2213 !/
2214 #ifdef W3_S
2215  INTEGER, SAVE :: IENT = 0
2216 #endif
2217 !/
2218 !/ ------------------------------------------------------------------- /
2219 !/
2220  integer ISEA, JSEA
2221  integer ierr
2222  real(rkind), intent(out) :: SXX_t(NSEA), SXY_t(NSEA), SYY_t(NSEA)
2223  real(rkind) :: SXX_p(NSEAL), SXY_p(NSEAL), SYY_p(NSEAL)
2224  real(rkind), allocatable :: rVect(:)
2225  integer IPROC, NSEAL_loc
2226 #ifdef W3_S
2227  CALL strace (ient, 'VA_SETUP_IOBPD')
2228 #endif
2229  DO isea=1,nseal
2230  sxx_p(isea)=sxx(isea)
2231  sxy_p(isea)=sxy(isea)
2232  syy_p(isea)=syy(isea)
2233  END DO
2234  IF (iaproc .eq. 1) THEN
2235  DO jsea=1,nseal
2236  isea=1 + (jsea-1)*naproc
2237  sxx_t(isea)=sxx_p(jsea)
2238  sxy_t(isea)=sxy_p(jsea)
2239  syy_t(isea)=syy_p(jsea)
2240  END DO
2241  DO iproc=2,naproc
2242  nseal_loc=1 + (nsea-iproc)/naproc
2243  allocate(rvect(nseal_loc))
2244  CALL mpi_recv(rvect,nseal_loc,rtype, iproc-1, 83, mpi_comm_wcmp, istatus, ierr)
2245  DO jsea=1,nseal_loc
2246  isea = iproc + (jsea-1)*naproc
2247  sxx_t(isea)=rvect(jsea)
2248  END DO
2249  CALL mpi_recv(rvect,nseal_loc,rtype, iproc-1, 89, mpi_comm_wcmp, istatus, ierr)
2250  DO jsea=1,nseal_loc
2251  isea = iproc + (jsea-1)*naproc
2252  sxy_t(isea)=rvect(jsea)
2253  END DO
2254  CALL mpi_recv(rvect,nseal_loc,rtype, iproc-1, 97, mpi_comm_wcmp, istatus, ierr)
2255  DO jsea=1,nseal_loc
2256  isea = iproc + (jsea-1)*naproc
2257  syy_t(isea)=rvect(jsea)
2258  END DO
2259  deallocate(rvect)
2260  END DO
2261  ELSE
2262  CALL mpi_send(sxx_p,nseal,rtype, 0, 83, mpi_comm_wcmp, ierr)
2263  CALL mpi_send(sxy_p,nseal,rtype, 0, 83, mpi_comm_wcmp, ierr)
2264  CALL mpi_send(syy_p,nseal,rtype, 0, 83, mpi_comm_wcmp, ierr)
2265  END IF

References w3odatmd::iaproc, yowdatapool::istatus, w3adatmd::mpi_comm_wcmp, w3odatmd::naproc, w3gdatmd::nsea, w3gdatmd::nseal, yowdatapool::rtype, w3servmd::strace(), w3adatmd::sxx, w3adatmd::sxy, and w3adatmd::syy.

Referenced by fd_wave_setup_computation().

◆ fd_compute_diff()

subroutine w3wavset::fd_compute_diff ( integer, intent(in)  IEDGE,
integer, intent(in)  ISEA,
real(rkind), intent(inout)  UGRAD,
real(rkind), intent(inout)  VGRAD,
real(rkind), intent(inout)  dist 
)

Differences on FD grids.

Parameters
[in]IEDGE
[in]ISEA
[in,out]UGRAD
[in,out]VGRAD
[in,out]dist
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2434 of file w3wavset.F90.

2434 !/
2435 !/ +-----------------------------------+
2436 !/ | WAVEWATCH III NOAA/NCEP |
2437 !/ | |
2438 !/ | Mathieu Dutour-Sikiric (IRB) |
2439 !/ | Aron Roland (BGS IT&E GmbH) |
2440 !/ | |
2441 !/ | FORTRAN 90 |
2442 !/ | Last update : 01-Mai-2018 |
2443 !/ +-----------------------------------+
2444 !/
2445 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2446 !/
2447 ! 1. Purpose : differences on FD grids
2448 ! 2. Method :
2449 ! 3. Parameters :
2450 !
2451 ! Parameter list
2452 ! ----------------------------------------------------------------
2453 ! ----------------------------------------------------------------
2454 !
2455 ! 4. Subroutines used :
2456 !
2457 ! Name Type Module Description
2458 ! ----------------------------------------------------------------
2459 ! STRACE Subr. W3SERVMD Subroutine tracing.
2460 ! ----------------------------------------------------------------
2461 !
2462 ! 5. Called by :
2463 !
2464 ! Name Type Module Description
2465 ! ----------------------------------------------------------------
2466 ! ----------------------------------------------------------------
2467 !
2468 ! 6. Error messages :
2469 ! 7. Remarks
2470 ! 8. Structure :
2471 ! 9. Switches :
2472 !
2473 ! !/S Enable subroutine tracing.
2474 !
2475 ! 10. Source code :
2476 !
2477 !/ ------------------------------------------------------------------- /
2478 #ifdef W3_S
2479  USE w3servmd, ONLY: strace
2480 #endif
2481 !
2482  USE w3gdatmd, ONLY: mapsf, edges
2483  USE w3gdatmd, ONLY: xgrd, ygrd
2484  IMPLICIT NONE
2485 !/
2486 !/ ------------------------------------------------------------------- /
2487 !/ Parameter list
2488 !/
2489 !/ ------------------------------------------------------------------- /
2490 !/ Local PARAMETERs
2491 !/
2492 #ifdef W3_S
2493  INTEGER, SAVE :: IENT = 0
2494 #endif
2495 !/
2496 !/ ------------------------------------------------------------------- /
2497 !/
2498  INTEGER, intent(in) :: IEDGE, ISEA
2499  REAL(rkind), intent(inout) :: UGRAD, VGRAD, dist
2500  REAL(rkind) :: h
2501  integer I2, I3, IP1, IP2, IP3
2502  integer IX1, IY1, IX2, IY2
2503  integer ISEA1, ISEA2
2504  REAL(rkind) deltaX, deltaY
2505 #ifdef W3_S
2506  CALL strace (ient, 'VA_SETUP_IOBPD')
2507 #endif
2508  !
2509  isea1=edges(iedge,1)
2510  isea2=edges(iedge,2)
2511  ix1=mapsf(isea1,1)
2512  iy1=mapsf(isea1,2)
2513  ix2=mapsf(isea2,1)
2514  iy2=mapsf(isea2,2)
2515  deltax=xgrd(ix1,iy1) - xgrd(ix2,iy2)
2516  deltay=ygrd(ix1,iy1) - ygrd(ix2,iy2)
2517  dist=sqrt(deltax*deltax + deltay*deltay)
2518  IF (isea .eq. isea1) THEN
2519  ugrad= deltax/dist
2520  vgrad= deltay/dist
2521  ELSE
2522  ugrad=-deltax/dist
2523  vgrad=-deltay/dist
2524  END IF

References w3gdatmd::edges, w3gdatmd::mapsf, w3servmd::strace(), w3gdatmd::xgrd, and w3gdatmd::ygrd.

Referenced by fd_wave_setup_compute_system().

◆ fd_compute_lh_stress()

subroutine w3wavset::fd_compute_lh_stress ( real(rkind), dimension(nsea), intent(in)  SXX_t,
real(rkind), dimension(nsea), intent(in)  SXY_t,
real(rkind), dimension(nsea), intent(in)  SYY_t,
real(rkind), dimension(nsea), intent(out)  FX,
real(rkind), dimension(nsea), intent(out)  FY 
)

Setup fluxes.

Parameters
[in]SXX_t
[in]SXY_t
[in]SYY_t
[out]FX
[out]FY
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2282 of file w3wavset.F90.

2282 !/
2283 !/ +-----------------------------------+
2284 !/ | WAVEWATCH III NOAA/NCEP |
2285 !/ | |
2286 !/ | Mathieu Dutour-Sikiric (IRB) |
2287 !/ | Aron Roland (BGS IT&E GmbH) |
2288 !/ | |
2289 !/ | FORTRAN 90 |
2290 !/ | Last update : 01-Mai-2018 |
2291 !/ +-----------------------------------+
2292 !/
2293 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2294 !/
2295 ! 1. Purpose : setup fluxes
2296 ! 2. Method :
2297 ! 3. Parameters :
2298 !
2299 ! Parameter list
2300 ! ----------------------------------------------------------------
2301 ! ----------------------------------------------------------------
2302 !
2303 ! 4. Subroutines used :
2304 !
2305 ! Name Type Module Description
2306 ! ----------------------------------------------------------------
2307 ! STRACE Subr. W3SERVMD Subroutine tracing.
2308 ! ----------------------------------------------------------------
2309 !
2310 ! 5. Called by :
2311 !
2312 ! Name Type Module Description
2313 ! ----------------------------------------------------------------
2314 ! ----------------------------------------------------------------
2315 !
2316 ! 6. Error messages :
2317 ! 7. Remarks
2318 ! 8. Structure :
2319 ! 9. Switches :
2320 !
2321 ! !/S Enable subroutine tracing.
2322 !
2323 ! 10. Source code :
2324 !
2325 !/ ------------------------------------------------------------------- /
2326 #ifdef W3_S
2327  USE w3servmd, ONLY: strace
2328 #endif
2329 !
2330  USE w3gdatmd, ONLY: nx, ny, nsea, neigh
2331  USE w3adatmd, ONLY: sxx, sxy, syy
2332  IMPLICIT NONE
2333 !/
2334 !/ ------------------------------------------------------------------- /
2335 !/ Parameter list
2336 !/
2337 !/ ------------------------------------------------------------------- /
2338 !/ Local PARAMETERs
2339 !/
2340 #ifdef W3_S
2341  INTEGER, SAVE :: IENT = 0
2342 #endif
2343 !/
2344 !/ ------------------------------------------------------------------- /
2345 !/
2346  real(rkind), intent(in) :: SXX_t(NSEA), SXY_t(NSEA), SYY_t(NSEA)
2347  real(rkind), intent(out) :: FX(NSEA), FY(NSEA)
2348  REAL(rkind) :: h
2349  REAL(rkind) :: SXX_X, SXX_Y
2350  REAL(rkind) :: SXY_X, SXY_Y
2351  REAL(rkind) :: SYY_X, SYY_Y
2352  REAL(rkind) :: eFX, eFY
2353  REAL(rkind) :: UGRAD, VGRAD
2354  INTEGER IE, I1, I2, I3, IP1, IP2, IP3
2355  integer ISEA, JSEA1, JSEA2, JSEA3, JSEA4
2356  integer NeighMat(4,2)
2357  real(rkind) dist_X, dist_Y
2358 #ifdef W3_S
2359  CALL strace (ient, 'VA_SETUP_IOBPD')
2360 #endif
2361  !
2362  neighmat(1,1)=1
2363  neighmat(1,2)=0
2364  neighmat(2,1)=-1
2365  neighmat(2,2)=0
2366  neighmat(3,1)=0
2367  neighmat(3,2)=1
2368  neighmat(4,1)=0
2369  neighmat(4,2)=-1
2370  fx=0
2371  fy=0
2372  DO isea=1,nsea
2373  jsea1=neigh(isea,1)
2374  jsea2=neigh(isea,2)
2375  jsea3=neigh(isea,3)
2376  jsea4=neigh(isea,4)
2377  sxx_x=0
2378  sxx_y=0
2379  sxy_x=0
2380  sxy_y=0
2381  syy_x=0
2382  syy_y=0
2383  IF ((jsea1 .gt. 0).and.(jsea2 .gt. 0)) THEN
2384  sxx_x=(sxx(jsea1) - sxx(jsea2))/(2*dist_x)
2385  sxy_x=(sxy(jsea1) - sxy(jsea2))/(2*dist_x)
2386  syy_x=(sxy(jsea1) - syy(jsea2))/(2*dist_x)
2387  END IF
2388  IF ((jsea1 .gt. 0).and.(jsea2 .eq. 0)) THEN
2389  sxx_x=(sxx(jsea1) - sxx(isea ))/dist_x
2390  sxy_x=(sxy(jsea1) - sxy(isea ))/dist_x
2391  syy_x=(sxy(jsea1) - syy(isea ))/dist_x
2392  END IF
2393  IF ((jsea1 .eq. 0).and.(jsea2 .gt. 0)) THEN
2394  sxx_x=(sxx(isea ) - sxx(jsea2))/dist_x
2395  sxy_x=(sxy(isea ) - sxy(jsea2))/dist_x
2396  syy_x=(sxy(isea ) - syy(jsea2))/dist_x
2397  END IF
2398  IF ((jsea3 .gt. 0).and.(jsea4 .gt. 0)) THEN
2399  sxx_x=(sxx(jsea3) - sxx(jsea4))/(2*dist_y)
2400  sxy_x=(sxy(jsea3) - sxy(jsea4))/(2*dist_y)
2401  syy_x=(sxy(jsea3) - syy(jsea4))/(2*dist_y)
2402  END IF
2403  IF ((jsea3 .eq. 0).and.(jsea4 .gt. 0)) THEN
2404  sxx_x=(sxx(isea ) - sxx(jsea4))/dist_y
2405  sxy_x=(sxy(isea ) - sxy(jsea4))/dist_y
2406  syy_x=(sxy(isea ) - syy(jsea4))/dist_y
2407  END IF
2408  IF ((jsea3 .gt. 0).and.(jsea4 .gt. 0)) THEN
2409  sxx_x=(sxx(jsea3) - sxx(isea ))/dist_y
2410  sxy_x=(sxy(jsea3) - sxy(isea ))/dist_y
2411  syy_x=(sxy(jsea3) - syy(isea ))/dist_y
2412  END IF
2413  efx=-sxx_x - sxy_y
2414  efy=-syy_y - sxy_x
2415  fx(isea)=efx
2416  fy(isea)=efy
2417  END DO

References w3gdatmd::neigh, w3gdatmd::nsea, w3gdatmd::nx, w3gdatmd::ny, w3servmd::strace(), w3adatmd::sxx, w3adatmd::sxy, and w3adatmd::syy.

Referenced by fd_wave_setup_computation().

◆ fd_set_meanvalue_to_zero()

subroutine w3wavset::fd_set_meanvalue_to_zero ( real(rkind), dimension(nx), intent(inout)  TheVar)

Set mean value.

Parameters
[in,out]TheVar
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2864 of file w3wavset.F90.

2864 !/
2865 !/ +-----------------------------------+
2866 !/ | WAVEWATCH III NOAA/NCEP |
2867 !/ | |
2868 !/ | Mathieu Dutour-Sikiric (IRB) |
2869 !/ | Aron Roland (BGS IT&E GmbH) |
2870 !/ | |
2871 !/ | FORTRAN 90 |
2872 !/ | Last update : 01-Mai-2018 |
2873 !/ +-----------------------------------+
2874 !/
2875 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2876 !/
2877 ! 1. Purpose : set meanvalue
2878 ! 2. Method :
2879 ! 3. Parameters :
2880 !
2881 ! Parameter list
2882 ! ----------------------------------------------------------------
2883 ! ----------------------------------------------------------------
2884 !
2885 ! 4. Subroutines used :
2886 !
2887 ! Name Type Module Description
2888 ! ----------------------------------------------------------------
2889 ! STRACE Subr. W3SERVMD Subroutine tracing.
2890 ! ----------------------------------------------------------------
2891 !
2892 ! 5. Called by :
2893 !
2894 ! Name Type Module Description
2895 ! ----------------------------------------------------------------
2896 ! ----------------------------------------------------------------
2897 !
2898 ! 6. Error messages :
2899 ! 7. Remarks
2900 ! 8. Structure :
2901 ! 9. Switches :
2902 !
2903 ! !/S Enable subroutine tracing.
2904 !
2905 ! 10. Source code :
2906 !
2907 !/ ------------------------------------------------------------------- /
2908 #ifdef W3_S
2909  USE w3servmd, ONLY: strace
2910 #endif
2911 !
2912  USE w3gdatmd, ONLY: nx, si
2913  IMPLICIT NONE
2914 !/
2915 !/ ------------------------------------------------------------------- /
2916 !/ Parameter list
2917 !/
2918 !/ ------------------------------------------------------------------- /
2919 !/ Local PARAMETERs
2920 !/
2921 #ifdef W3_S
2922  INTEGER, SAVE :: IENT = 0
2923 #endif
2924 !/
2925 !/ ------------------------------------------------------------------- /
2926 !/
2927  real(rkind), intent(inout) :: TheVar(NX)
2928  real(rkind) :: SUM_SI_Var, SUM_SI, TheMean
2929  INTEGER IP
2930 #ifdef W3_S
2931  CALL strace (ient, 'VA_SETUP_IOBPD')
2932 #endif
2933  sum_si_var=0
2934  sum_si=0
2935  DO ip=1,nx
2936  sum_si_var = sum_si_var + si(ip)*thevar(ip)
2937  sum_si = sum_si + si(ip)
2938  END DO
2939  themean=sum_si_var/sum_si
2940  DO ip=1,nx
2941  thevar(ip)=thevar(ip) - themean
2942  END DO

References w3gdatmd::nx, w3gdatmd::si, and w3servmd::strace().

Referenced by fd_wave_setup_computation().

◆ fd_wave_setup_apply_fct()

subroutine w3wavset::fd_wave_setup_apply_fct ( real(rkind), dimension(nnz), intent(in)  ASPAR,
real(rkind), dimension(nsea), intent(in)  TheIn,
real(rkind), dimension(nsea), intent(out)  TheOut 
)

Compute off diagonal for FD grids.

Parameters
[in]ASPAR
[in]TheIn
[out]TheOut
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1946 of file w3wavset.F90.

1946 !/
1947 !/ +-----------------------------------+
1948 !/ | WAVEWATCH III NOAA/NCEP |
1949 !/ | |
1950 !/ | Mathieu Dutour-Sikiric (IRB) |
1951 !/ | Aron Roland (BGS IT&E GmbH) |
1952 !/ | |
1953 !/ | FORTRAN 90 |
1954 !/ | Last update : 01-Mai-2018 |
1955 !/ +-----------------------------------+
1956 !/
1957 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1958 !/
1959 ! 1. Purpose : comp. off diagonal for FD grids
1960 ! 2. Method :
1961 ! 3. Parameters :
1962 !
1963 ! Parameter list
1964 ! ----------------------------------------------------------------
1965 ! ----------------------------------------------------------------
1966 !
1967 ! 4. Subroutines used :
1968 !
1969 ! Name Type Module Description
1970 ! ----------------------------------------------------------------
1971 ! STRACE Subr. W3SERVMD Subroutine tracing.
1972 ! ----------------------------------------------------------------
1973 !
1974 ! 5. Called by :
1975 !
1976 ! Name Type Module Description
1977 ! ----------------------------------------------------------------
1978 ! ----------------------------------------------------------------
1979 !
1980 ! 6. Error messages :
1981 ! 7. Remarks
1982 ! 8. Structure :
1983 ! 9. Switches :
1984 !
1985 ! !/S Enable subroutine tracing.
1986 !
1987 ! 10. Source code :
1988 !
1989 !/ ------------------------------------------------------------------- /
1990 #ifdef W3_S
1991  USE w3servmd, ONLY: strace
1992 #endif
1993 !
1994  USE w3gdatmd, ONLY: nx, nnz, iaa, jaa, nsea
1995  use yownodepool, only: pdlib_ia, pdlib_ja
1996  IMPLICIT NONE
1997 !/
1998 !/ ------------------------------------------------------------------- /
1999 !/ Parameter list
2000 !/
2001 !/ ------------------------------------------------------------------- /
2002 !/ Local PARAMETERs
2003 !/
2004 #ifdef W3_S
2005  INTEGER, SAVE :: IENT = 0
2006 #endif
2007 !/
2008 !/ ------------------------------------------------------------------- /
2009 !/
2010  REAL(rkind), intent(in) :: ASPAR(NNZ)
2011  REAL(rkind), intent(in) :: TheIn(NSEA)
2012  REAL(rkind), intent(out) :: TheOut(NSEA)
2013  integer IP, J, JP
2014  REAL(rkind) :: eCoeff
2015 #ifdef W3_S
2016  CALL strace (ient, 'VA_SETUP_IOBPD')
2017 #endif
2018  theout=0
2019  DO ip=1,nsea
2020  DO j=pdlib_ia(ip),pdlib_ia(ip+1)-1
2021  jp=pdlib_ja(j)
2022  ecoeff=aspar(j)
2023  theout(ip)=theout(ip) + ecoeff*thein(jp)
2024  END DO
2025  END DO

References w3gdatmd::iaa, w3gdatmd::jaa, w3gdatmd::nnz, w3gdatmd::nsea, w3gdatmd::nx, yownodepool::pdlib_ia, yownodepool::pdlib_ja, and w3servmd::strace().

Referenced by fd_wave_setup_solve_poisson_neumann_dir().

◆ fd_wave_setup_apply_precond()

subroutine w3wavset::fd_wave_setup_apply_precond ( real(rkind), dimension(pdlib_nnz), intent(in)  ASPAR,
real(rkind), dimension(nsea), intent(in)  TheIn,
real(rkind), dimension(nsea), intent(out)  TheOut 
)

Preconditioning for FD grids.

Parameters
[in]ASPAR
[in]TheIn
[out]TheOut
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2040 of file w3wavset.F90.

2040 !/
2041 !/ +-----------------------------------+
2042 !/ | WAVEWATCH III NOAA/NCEP |
2043 !/ | |
2044 !/ | Mathieu Dutour-Sikiric (IRB) |
2045 !/ | Aron Roland (BGS IT&E GmbH) |
2046 !/ | |
2047 !/ | FORTRAN 90 |
2048 !/ | Last update : 01-Mai-2018 |
2049 !/ +-----------------------------------+
2050 !/
2051 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2052 !/
2053 ! 1. Purpose : Precond. for FD grids
2054 ! 2. Method :
2055 ! 3. Parameters :
2056 !
2057 ! Parameter list
2058 ! ----------------------------------------------------------------
2059 ! ----------------------------------------------------------------
2060 !
2061 ! 4. Subroutines used :
2062 !
2063 ! Name Type Module Description
2064 ! ----------------------------------------------------------------
2065 ! STRACE Subr. W3SERVMD Subroutine tracing.
2066 ! ----------------------------------------------------------------
2067 !
2068 ! 5. Called by :
2069 !
2070 ! Name Type Module Description
2071 ! ----------------------------------------------------------------
2072 ! ----------------------------------------------------------------
2073 !
2074 ! 6. Error messages :
2075 ! 7. Remarks
2076 ! 8. Structure :
2077 ! 9. Switches :
2078 !
2079 ! !/S Enable subroutine tracing.
2080 !
2081 ! 10. Source code :
2082 !
2083 !/ ------------------------------------------------------------------- /
2084 #ifdef W3_S
2085  USE w3servmd, ONLY: strace
2086 #endif
2087 !
2089  USE w3gdatmd, ONLY: nsea
2090  IMPLICIT NONE
2091 !/
2092 !/ ------------------------------------------------------------------- /
2093 !/ Parameter list
2094 !/
2095 !/ ------------------------------------------------------------------- /
2096 !/ Local PARAMETERs
2097 !/
2098 #ifdef W3_S
2099  INTEGER, SAVE :: IENT = 0
2100 #endif
2101 !/
2102 !/ ------------------------------------------------------------------- /
2103 !/
2104  REAL(rkind), intent(in) :: ASPAR(PDLIB_NNZ)
2105  REAL(rkind), intent(in) :: TheIn(NSEA)
2106  REAL(rkind), intent(out) :: TheOut(NSEA)
2107  integer IP, J1, J, JP, J2
2108  REAL(rkind) :: eCoeff
2109  INTEGER :: ThePrecond = 0
2110 #ifdef W3_S
2111  CALL strace (ient, 'VA_SETUP_IOBPD')
2112 #endif
2113  IF (theprecond .eq. 0) THEN
2114  theout=thein
2115  END IF
2116  IF (theprecond .eq. 1) THEN
2117  theout=0
2118  DO ip=1,nsea
2119  j1=pdlib_i_diag(ip)
2120  DO j=pdlib_ia(ip),pdlib_ia(ip+1)-1
2121  jp=pdlib_ja(j)
2122  IF (j .eq. j1) THEN
2123  ecoeff=1.0/aspar(j)
2124  ELSE
2125  j2=pdlib_i_diag(jp)
2126  ecoeff=-aspar(j) /(aspar(j1)*aspar(j2))
2127  END IF
2128  theout(ip)=theout(ip) + ecoeff*thein(jp)
2129  END DO
2130  END DO
2131  END IF
2132  IF (theprecond .eq. 2) THEN
2133 
2134  DO ip=1,nsea
2135  j=pdlib_i_diag(ip)
2136  theout(ip)=thein(ip)/aspar(j)
2137  END DO
2138  END IF

References w3gdatmd::nsea, yownodepool::pdlib_i_diag, yownodepool::pdlib_ia, yownodepool::pdlib_ja, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by fd_wave_setup_solve_poisson_neumann_dir().

◆ fd_wave_setup_computation()

subroutine w3wavset::fd_wave_setup_computation

Wave setup comp on FD grids.

Parameters
[in,out]TheVar
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2955 of file w3wavset.F90.

2955 !/
2956 !/ +-----------------------------------+
2957 !/ | WAVEWATCH III NOAA/NCEP |
2958 !/ | |
2959 !/ | Mathieu Dutour-Sikiric (IRB) |
2960 !/ | Aron Roland (BGS IT&E GmbH) |
2961 !/ | |
2962 !/ | FORTRAN 90 |
2963 !/ | Last update : 01-Mai-2018 |
2964 !/ +-----------------------------------+
2965 !/
2966 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2967 !/
2968 ! 1. Purpose : wave setup comp. on fd grids
2969 ! 2. Method :
2970 ! 3. Parameters :
2971 !
2972 ! Parameter list
2973 ! ----------------------------------------------------------------
2974 ! ----------------------------------------------------------------
2975 !
2976 ! 4. Subroutines used :
2977 !
2978 ! Name Type Module Description
2979 ! ----------------------------------------------------------------
2980 ! STRACE Subr. W3SERVMD Subroutine tracing.
2981 ! ----------------------------------------------------------------
2982 !
2983 ! 5. Called by :
2984 !
2985 ! Name Type Module Description
2986 ! ----------------------------------------------------------------
2987 ! ----------------------------------------------------------------
2988 !
2989 ! 6. Error messages :
2990 ! 7. Remarks
2991 ! 8. Structure :
2992 ! 9. Switches :
2993 !
2994 ! !/S Enable subroutine tracing.
2995 !
2996 ! 10. Source code :
2997 !
2998 !/ ------------------------------------------------------------------- /
2999 #ifdef W3_S
3000  USE w3servmd, ONLY: strace
3001 #endif
3002 !
3003  USE yownodepool, only: pdlib_nnz
3004  USE w3gdatmd, ONLY: nx, nsea, nseal
3005  USE w3wdatmd, ONLY: zeta_setup
3006  use yowdatapool, only: rtype, istatus
3007  USE w3adatmd, ONLY: mpi_comm_wcmp
3008  USE w3odatmd, only : iaproc, naproc
3009  IMPLICIT NONE
3010 !/
3011 !/ ------------------------------------------------------------------- /
3012 !/ Parameter list
3013 !/
3014 !/ ------------------------------------------------------------------- /
3015 !/ Local PARAMETERs
3016 !/
3017 #ifdef W3_S
3018  INTEGER, SAVE :: IENT = 0
3019 #endif
3020 !/
3021 !/ ------------------------------------------------------------------- /
3022 !/
3023  REAL(rkind) :: ZETA_WORK(NSEA)
3024  REAL(rkind) :: F_X(NSEA), F_Y(NSEA)
3025  REAL(rkind) :: ASPAR(PDLIB_NNZ), B(NX)
3026  INTEGER ISEA, IPROC
3027  real(rkind) :: SXX_t(NSEA), SXY_t(NSEA), SYY_t(NSEA)
3028  integer ierr
3029 #ifdef W3_DEBUGSTP
3030  real(rkind) max_val, min_val
3031 #endif
3032 #ifdef W3_S
3033  CALL strace (ient, 'VA_SETUP_IOBPD')
3034 #endif
3035  CALL fd_collect_sxx_xy_yy(sxx_t, sxy_t, syy_t)
3036  IF (iaproc .eq. 1) THEN
3037  CALL fd_compute_lh_stress(sxx_t, sxy_t, syy_t, f_x, f_y)
3038  DO isea=1,nsea
3039  zeta_work(isea)=zeta_setup(isea)
3040  END DO
3041  CALL fd_wave_setup_compute_system(aspar, b, f_x, f_y)
3042  CALL fd_wave_setup_solve_poisson_neumann_dir(aspar, b, zeta_work)
3043  CALL fd_set_meanvalue_to_zero(zeta_work)
3044  DO iproc=2,naproc
3045  CALL mpi_send(zeta_work,nsea,rtype, iproc-1, 23, mpi_comm_wcmp, ierr)
3046  END DO
3047  ELSE
3048  CALL mpi_recv(zeta_work,nseal,rtype, 0, 23, mpi_comm_wcmp, istatus, ierr)
3049  END IF
3050 #ifdef W3_DEBUGSTP
3051  max_val = zeta_work(isea)
3052  min_val = zeta_work(isea)
3053 #endif
3054  DO isea=1,nsea
3055  zeta_setup(isea)=zeta_work(isea)
3056 #ifdef W3_DEBUGSTP
3057  max_val = max(max_val, zeta_work(isea))
3058  min_val = min(min_val, zeta_work(isea))
3059 #endif
3060  END DO
3061 #ifdef W3_DEBUGSTP
3062  WRITE(740+iaproc,*) 'FD_WAVE_SETUP_COMPUTATION, max/min=', max_val, min_val
3063  FLUSH(740+iaproc)
3064 #endif

References fd_collect_sxx_xy_yy(), fd_compute_lh_stress(), fd_set_meanvalue_to_zero(), fd_wave_setup_compute_system(), fd_wave_setup_solve_poisson_neumann_dir(), w3odatmd::iaproc, yowdatapool::istatus, w3adatmd::mpi_comm_wcmp, w3odatmd::naproc, w3gdatmd::nsea, w3gdatmd::nseal, w3gdatmd::nx, yownodepool::pdlib_nnz, yowdatapool::rtype, w3servmd::strace(), and w3wdatmd::zeta_setup.

Referenced by wave_setup_computation().

◆ fd_wave_setup_compute_system()

subroutine w3wavset::fd_wave_setup_compute_system ( real(rkind), dimension(pdlib_nnz), intent(out)  ASPAR,
real(rkind), dimension(nx), intent(out)  B,
real(rkind), dimension(nsea), intent(in)  FX,
real(rkind), dimension(nsea), intent(in)  FY 
)

Setup matrix on FD grids.

Parameters
[out]ASPAR
[out]B
[in]FX
[in]FY
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2540 of file w3wavset.F90.

2540 !/
2541 !/ +-----------------------------------+
2542 !/ | WAVEWATCH III NOAA/NCEP |
2543 !/ | |
2544 !/ | Mathieu Dutour-Sikiric (IRB) |
2545 !/ | Aron Roland (BGS IT&E GmbH) |
2546 !/ | |
2547 !/ | FORTRAN 90 |
2548 !/ | Last update : 01-Mai-2018 |
2549 !/ +-----------------------------------+
2550 !/
2551 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2552 !/
2553 ! 1. Purpose : Setup matrix on FD grids
2554 ! 2. Method :
2555 ! 3. Parameters :
2556 !
2557 ! Parameter list
2558 ! ----------------------------------------------------------------
2559 ! ----------------------------------------------------------------
2560 !
2561 ! 4. Subroutines used :
2562 !
2563 ! Name Type Module Description
2564 ! ----------------------------------------------------------------
2565 ! STRACE Subr. W3SERVMD Subroutine tracing.
2566 ! ----------------------------------------------------------------
2567 !
2568 ! 5. Called by :
2569 !
2570 ! Name Type Module Description
2571 ! ----------------------------------------------------------------
2572 ! ----------------------------------------------------------------
2573 !
2574 ! 6. Error messages :
2575 ! 7. Remarks
2576 ! 8. Structure :
2577 ! 9. Switches :
2578 !
2579 ! !/S Enable subroutine tracing.
2580 !
2581 ! 10. Source code :
2582 !
2583 !/ ------------------------------------------------------------------- /
2584 #ifdef W3_S
2585  USE w3servmd, ONLY: strace
2586 #endif
2587 !
2588  USE yownodepool, only: pdlib_nnz
2589  USE w3gdatmd, ONLY: nx, ny, nsea, nbedge, edges
2590  USE w3adatmd, ONLY: dw
2591  IMPLICIT NONE
2592 !/
2593 !/ ------------------------------------------------------------------- /
2594 !/ Parameter list
2595 !/
2596 !/ ------------------------------------------------------------------- /
2597 !/ Local PARAMETERs
2598 !/
2599 #ifdef W3_S
2600  INTEGER, SAVE :: IENT = 0
2601 #endif
2602 !/
2603 !/ ------------------------------------------------------------------- /
2604 !/
2605  real(rkind), intent(in) :: FX(NSEA), FY(NSEA)
2606  real(rkind), intent(out) :: ASPAR(PDLIB_NNZ)
2607  real(rkind), intent(out) :: B(NX)
2608  INTEGER :: POS_TRICK(3,2), POS_SHIFT(3,3)
2609  integer I1, I2, I3, IP1, IP2, IP3
2610  integer IDX, IDX1, IDX2, IDX3
2611  INTEGER IE, IP, I, J, K, IPp, JPp
2612  real(rkind) :: eDep, eFX, eFY, eScal, eFact, eLen
2613  real(rkind) :: UGRAD, VGRAD, UGRAD1, VGRAD1, dist1, dist2
2614  INTEGER LIDX(2), KIDX(2), jdx
2615  INTEGER ISEAREL, JSEAREL, ISEA, JSEA, IEDGE
2616 #ifdef W3_S
2617  CALL strace (ient, 'VA_SETUP_IOBPD')
2618 #endif
2619  !
2620  aspar=0
2621  b=0
2622  DO iedge=1,nbedge
2623  isea=edges(iedge,1)
2624  jsea=edges(iedge,2)
2625  edep=(dw(isea) + dw(jsea))/2.0
2626  efx =(fx(isea) + fx(jsea))/2.0
2627  efy =(fy(isea) + fy(jsea))/2.0
2628  DO i=1,2
2629  isearel=edges(iedge,i)
2630  CALL fd_compute_diff(iedge, isearel, ugrad1, vgrad1, dist1)
2631  escal=ugrad1*efx + vgrad1*efy
2632  b(isearel) = b(isearel) + escal*dist1
2633  !
2634  DO j=1,2
2635  jsearel=edges(iedge,j)
2636  CALL fd_compute_diff(iedge, jsearel, ugrad, vgrad, dist2)
2637  escal=ugrad*ugrad1 + vgrad*vgrad1
2638  aspar(j)=aspar(j)+efact*escal
2639  END DO
2640  END DO
2641  END DO

References w3adatmd::dw, w3gdatmd::edges, fd_compute_diff(), w3gdatmd::nbedge, w3gdatmd::nsea, w3gdatmd::nx, w3gdatmd::ny, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by fd_wave_setup_computation().

◆ fd_wave_setup_scalar_prod()

subroutine w3wavset::fd_wave_setup_scalar_prod ( real(rkind), dimension(nx), intent(in)  V1,
real(rkind), dimension(nx), intent(in)  V2,
real(rkind), intent(inout)  eScal 
)

Scalar product.

Parameters
[in]V1
[in]V2
[in,out]eScal
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2656 of file w3wavset.F90.

2656 !/
2657 !/ +-----------------------------------+
2658 !/ | WAVEWATCH III NOAA/NCEP |
2659 !/ | |
2660 !/ | Mathieu Dutour-Sikiric (IRB) |
2661 !/ | Aron Roland (BGS IT&E GmbH) |
2662 !/ | |
2663 !/ | FORTRAN 90 |
2664 !/ | Last update : 01-Mai-2018 |
2665 !/ +-----------------------------------+
2666 !/
2667 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2668 !/
2669 ! 1. Purpose : scalar prod.
2670 ! 2. Method :
2671 ! 3. Parameters :
2672 !
2673 ! Parameter list
2674 ! ----------------------------------------------------------------
2675 ! ----------------------------------------------------------------
2676 !
2677 ! 4. Subroutines used :
2678 !
2679 ! Name Type Module Description
2680 ! ----------------------------------------------------------------
2681 ! STRACE Subr. W3SERVMD Subroutine tracing.
2682 ! ----------------------------------------------------------------
2683 !
2684 ! 5. Called by :
2685 !
2686 ! Name Type Module Description
2687 ! ----------------------------------------------------------------
2688 ! ----------------------------------------------------------------
2689 !
2690 ! 6. Error messages :
2691 ! 7. Remarks
2692 ! 8. Structure :
2693 ! 9. Switches :
2694 !
2695 ! !/S Enable subroutine tracing.
2696 !
2697 ! 10. Source code :
2698 !
2699 !/ ------------------------------------------------------------------- /
2700 #ifdef W3_S
2701  USE w3servmd, ONLY: strace
2702 #endif
2703 !
2704  USE w3gdatmd, ONLY: nx
2705  IMPLICIT NONE
2706 !/
2707 !/ ------------------------------------------------------------------- /
2708 !/ Parameter list
2709 !/
2710 !/ ------------------------------------------------------------------- /
2711 !/ Local PARAMETERs
2712 !/
2713 #ifdef W3_S
2714  INTEGER, SAVE :: IENT = 0
2715 #endif
2716 !/
2717 !/ ------------------------------------------------------------------- /
2718 !/
2719  real(rkind), intent(in) :: V1(NX), V2(NX)
2720  real(rkind), intent(inout) :: eScal
2721  integer IP
2722 #ifdef W3_S
2723  CALL strace (ient, 'VA_SETUP_IOBPD')
2724 #endif
2725  escal=0
2726  DO ip=1,nx
2727  escal=escal + v1(ip)*v2(ip)
2728  END DO

References w3gdatmd::nx, and w3servmd::strace().

Referenced by fd_wave_setup_solve_poisson_neumann_dir().

◆ fd_wave_setup_solve_poisson_neumann_dir()

subroutine w3wavset::fd_wave_setup_solve_poisson_neumann_dir ( real(rkind), dimension(pdlib_nnz), intent(in)  ASPAR,
real(rkind), dimension(nx), intent(in)  B,
real(rkind), dimension(nx), intent(out)  TheOut 
)

Poisson solver on FD grids.

Parameters
[in]ASPAR
[in]B
[out]TheOut
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 2743 of file w3wavset.F90.

2743 !/
2744 !/ +-----------------------------------+
2745 !/ | WAVEWATCH III NOAA/NCEP |
2746 !/ | |
2747 !/ | Mathieu Dutour-Sikiric (IRB) |
2748 !/ | Aron Roland (BGS IT&E GmbH) |
2749 !/ | |
2750 !/ | FORTRAN 90 |
2751 !/ | Last update : 01-Mai-2018 |
2752 !/ +-----------------------------------+
2753 !/
2754 !/ 01-Mai-2018 : Origination. ( version 6.04 )
2755 !/
2756 ! 1. Purpose : possoin solver on fd grids
2757 ! 2. Method :
2758 ! 3. Parameters :
2759 !
2760 ! Parameter list
2761 ! ----------------------------------------------------------------
2762 ! ----------------------------------------------------------------
2763 !
2764 ! 4. Subroutines used :
2765 !
2766 ! Name Type Module Description
2767 ! ----------------------------------------------------------------
2768 ! STRACE Subr. W3SERVMD Subroutine tracing.
2769 ! ----------------------------------------------------------------
2770 !
2771 ! 5. Called by :
2772 !
2773 ! Name Type Module Description
2774 ! ----------------------------------------------------------------
2775 ! ----------------------------------------------------------------
2776 !
2777 ! 6. Error messages :
2778 ! 7. Remarks
2779 ! 8. Structure :
2780 ! 9. Switches :
2781 !
2782 ! !/S Enable subroutine tracing.
2783 !
2784 ! 10. Source code :
2785 !
2786 !/ ------------------------------------------------------------------- /
2787 #ifdef W3_S
2788  USE w3servmd, ONLY: strace
2789 #endif
2790 !
2791  USE yownodepool, only: pdlib_nnz
2792  USE w3gdatmd, ONLY: nx
2793  IMPLICIT NONE
2794 !/
2795 !/ ------------------------------------------------------------------- /
2796 !/ Parameter list
2797 !/
2798 !/ ------------------------------------------------------------------- /
2799 !/ Local PARAMETERs
2800 !/
2801 #ifdef W3_S
2802  INTEGER, SAVE :: IENT = 0
2803 #endif
2804 !/
2805 !/ ------------------------------------------------------------------- /
2806 !/
2807  real(rkind), intent(in) :: ASPAR(PDLIB_NNZ)
2808  real(rkind), intent(in) :: B(NX)
2809  real(rkind), intent(out) :: TheOut(NX)
2810  real(rkind) :: V_X(NX), V_R(NX), V_Z(NX), V_P(NX), V_Y(NX)
2811  real(rkind) :: uO, uN, alphaV, h1, h2
2812  real(rkind) :: eNorm, beta
2813  real(rkind) :: SOLVERTHR
2814  integer IP, nbIter
2815 #ifdef W3_S
2816  CALL strace (ient, 'VA_SETUP_IOBPD')
2817 #endif
2818  solverthr=0.00000001
2819  nbiter=0
2820  v_x=0
2821  v_r=b
2822  CALL fd_wave_setup_apply_precond(aspar, v_r, v_z)
2823  v_p=v_z
2824  CALL fd_wave_setup_scalar_prod(v_z, v_r, uo)
2825  DO
2826  nbiter=nbiter + 1
2827  CALL fd_wave_setup_apply_fct(aspar, v_p, v_y)
2828  CALL fd_wave_setup_scalar_prod(v_p, v_y, h2)
2829  alphav=uo/h2
2830  !
2831  DO ip=1,nx
2832  v_x(ip) = v_x(ip) + alphav * v_p(ip)
2833  v_r(ip) = v_r(ip) - alphav * v_y(ip)
2834  END DO
2835  !
2836  CALL fd_wave_setup_scalar_prod(v_r, v_r, enorm)
2837  IF (enorm .le. solverthr) THEN
2838  EXIT
2839  END IF
2840  !
2841  CALL fd_wave_setup_apply_precond(aspar, v_r, v_z)
2842  CALL fd_wave_setup_scalar_prod(v_z, v_r, un)
2843  !
2844  beta=un/uo
2845  uo=un
2846  !
2847  DO ip=1,nx
2848  v_p(ip)=v_z(ip) + beta * v_p(ip)
2849  END DO
2850  END DO
2851  theout=v_x

References fd_wave_setup_apply_fct(), fd_wave_setup_apply_precond(), fd_wave_setup_scalar_prod(), w3gdatmd::nx, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by fd_wave_setup_computation().

◆ preparation_fd_scheme()

subroutine w3wavset::preparation_fd_scheme ( integer, intent(in)  IMOD)

Wave setup for FD grids.

Parameters
[in]IMOD
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1778 of file w3wavset.F90.

1778 !/
1779 !/ +-----------------------------------+
1780 !/ | WAVEWATCH III NOAA/NCEP |
1781 !/ | |
1782 !/ | Mathieu Dutour-Sikiric (IRB) |
1783 !/ | Aron Roland (BGS IT&E GmbH) |
1784 !/ | |
1785 !/ | FORTRAN 90 |
1786 !/ | Last update : 01-Mai-2018 |
1787 !/ +-----------------------------------+
1788 !/
1789 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1790 !/
1791 ! 1. Purpose : Wave setup for FD grids
1792 ! 2. Method :
1793 ! 3. Parameters :
1794 !
1795 ! Parameter list
1796 ! ----------------------------------------------------------------
1797 ! ----------------------------------------------------------------
1798 !
1799 ! 4. Subroutines used :
1800 !
1801 ! Name Type Module Description
1802 ! ----------------------------------------------------------------
1803 ! STRACE Subr. W3SERVMD Subroutine tracing.
1804 ! ----------------------------------------------------------------
1805 !
1806 ! 5. Called by :
1807 !
1808 ! Name Type Module Description
1809 ! ----------------------------------------------------------------
1810 ! ----------------------------------------------------------------
1811 !
1812 ! 6. Error messages :
1813 ! 7. Remarks
1814 ! 8. Structure :
1815 ! 9. Switches :
1816 !
1817 ! !/S Enable subroutine tracing.
1818 !
1819 ! 10. Source code :
1820 !
1821 !/ ------------------------------------------------------------------- /
1822 #ifdef W3_S
1823  USE w3servmd, ONLY: strace
1824 #endif
1825 !
1827  USE w3gdatmd, ONLY: nx, ny, nsea, mapsf, grids
1828  IMPLICIT NONE
1829 !/
1830 !/ ------------------------------------------------------------------- /
1831 !/ Parameter list
1832 !/
1833 !/ ------------------------------------------------------------------- /
1834 !/ Local PARAMETERs
1835 !/
1836 #ifdef W3_S
1837  INTEGER, SAVE :: IENT = 0
1838 #endif
1839 !/
1840 !/ ------------------------------------------------------------------- /
1841 !/
1842  integer, intent(in) :: IMOD
1843  integer IN, ISEA, nbEdge
1844  integer IX, IY, idx
1845  integer NeighMat(4,2)
1846  integer, allocatable :: STAT_SeaLand(:,:)
1847  integer, allocatable :: EDGES(:,:)
1848  integer IXN, JXN, JSEA, J
1849 #ifdef W3_S
1850  CALL strace (ient, 'VA_SETUP_IOBPD')
1851 #endif
1852  !
1853  allocate(grids(imod)%NEIGH(nsea,4))
1854  grids(imod)%NEIGH=0
1855  allocate(stat_sealand(nx,ny))
1856  stat_sealand=0
1857  DO isea=1,nsea
1858  ix=mapsf(isea,1)
1859  iy=mapsf(isea,2)
1860  stat_sealand(ix,iy)=isea
1861  END DO
1862  neighmat(1,1)=1
1863  neighmat(1,2)=0
1864  neighmat(2,1)=-1
1865  neighmat(2,2)=0
1866  neighmat(3,1)=0
1867  neighmat(3,2)=1
1868  neighmat(4,1)=0
1869  neighmat(4,2)=-1
1870  nbedge=0
1871  pdlib_nnz=0
1872  DO isea=1,nsea
1873  ix=mapsf(isea,1)
1874  iy=mapsf(isea,2)
1875  idx=0
1876  DO in=1,4
1877  ixn=ix+neighmat(in,1)
1878  jxn=ix+neighmat(in,2)
1879  jsea=stat_sealand(ixn,jxn)
1880  IF (jsea .gt. 0) THEN
1881  idx=idx+1
1882  grids(imod)%NEIGH(isea,idx)=jsea
1883  IF (jsea < isea) THEN
1884  nbedge=nbedge+1
1885  END IF
1887  END IF
1888  END DO
1890  END DO
1891  !
1892  grids(imod)%NBEDGE=nbedge
1893  ALLOCATE(grids(imod)%EDGES(nbedge,2))
1894  idx=0
1895  DO isea=1,nsea
1896  ix=mapsf(isea,1)
1897  iy=mapsf(isea,2)
1898  DO in=1,4
1899  ixn=ix+neighmat(in,1)
1900  jxn=ix+neighmat(in,2)
1901  jsea=stat_sealand(ixn,jxn)
1902  IF (jsea .gt. 0) THEN
1903  IF (jsea < isea) THEN
1904  idx=idx+1
1905  grids(imod)%EDGES(idx,1)=jsea
1906  grids(imod)%EDGES(idx,2)=isea
1907  END IF
1908  END IF
1909  END DO
1910  END DO
1911  !
1912  ALLOCATE(pdlib_ia(nsea+1))
1913  ALLOCATE(pdlib_ja(pdlib_nnz))
1914  ALLOCATE(pdlib_i_diag(nsea))
1915  pdlib_ia(1)=1
1916  j=0
1917  DO isea=1,nsea
1918  DO in=1,4
1919  ixn=ix+neighmat(in,1)
1920  jxn=ix+neighmat(in,2)
1921  jsea=stat_sealand(ixn,jxn)
1922  IF (jsea .gt. 0) THEN
1923  j=j+1
1924  pdlib_ja(j)=jsea
1925  END IF
1926  END DO
1927  j=j+1
1928  pdlib_ja(j)=isea
1929  pdlib_i_diag(isea)=j
1930  pdlib_ia(isea+1)=j+1
1931  END DO

References w3gdatmd::grids, w3gdatmd::mapsf, w3gdatmd::nsea, w3gdatmd::nx, w3gdatmd::ny, yownodepool::pdlib_i_diag, yownodepool::pdlib_ia, yownodepool::pdlib_ja, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by w3initmd::w3init().

◆ trig_compute_diff()

subroutine w3wavset::trig_compute_diff ( integer, intent(in)  IE,
integer, intent(in)  I1,
real(rkind), intent(inout)  UGRAD,
real(rkind), intent(inout)  VGRAD 
)

Differentiate other way around.

Parameters
[in]IE
[in]I1
[in,out]UGRAD
[in,out]VGRAD
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 599 of file w3wavset.F90.

599 !/
600 !/ +-----------------------------------+
601 !/ | WAVEWATCH III NOAA/NCEP |
602 !/ | |
603 !/ | Mathieu Dutour-Sikiric (IRB) |
604 !/ | Aron Roland (BGS IT&E GmbH) |
605 !/ | |
606 !/ | FORTRAN 90 |
607 !/ | Last update : 01-Mai-2018 |
608 !/ +-----------------------------------+
609 !/
610 !/ 01-Mai-2018 : Origination. ( version 6.04 )
611 !/
612 ! 1. Purpose : differentiate other way around ...
613 ! 2. Method :
614 ! 3. Parameters :
615 !
616 ! Parameter list
617 ! ----------------------------------------------------------------
618 ! ----------------------------------------------------------------
619 !
620 ! 4. Subroutines used :
621 !
622 ! Name Type Module Description
623 ! ----------------------------------------------------------------
624 ! STRACE Subr. W3SERVMD Subroutine tracing.
625 ! ----------------------------------------------------------------
626 !
627 ! 5. Called by :
628 !
629 ! Name Type Module Description
630 ! ----------------------------------------------------------------
631 ! ----------------------------------------------------------------
632 !
633 ! 6. Error messages :
634 ! 7. Remarks
635 ! 8. Structure :
636 ! 9. Switches :
637 !
638 ! !/S Enable subroutine tracing.
639 !
640 ! 10. Source code :
641 !
642 !/ ------------------------------------------------------------------- /
643 #ifdef W3_S
644  USE w3servmd, ONLY: strace
645 #endif
646 !
647  use yowelementpool, only: ine
648  use yownodepool, only: x, y, pdlib_tria
649  IMPLICIT NONE
650 !/
651 !/ ------------------------------------------------------------------- /
652 !/ Parameter list
653 !/
654 !/ ------------------------------------------------------------------- /
655 !/ Local PARAMETERs
656 !/
657 #ifdef W3_S
658  INTEGER, SAVE :: IENT = 0
659 #endif
660 !/
661 !/ ------------------------------------------------------------------- /
662 !/
663  INTEGER, intent(in) :: IE, I1
664  REAL(rkind), intent(inout) :: UGRAD, VGRAD
665  REAL(rkind) :: h
666  integer I2, I3, IP1, IP2, IP3
667  INTEGER :: POS_TRICK(3,2)
668 #ifdef W3_S
669  CALL strace (ient, 'VA_SETUP_IOBPD')
670 #endif
671  pos_trick(1,1) = 2
672  pos_trick(1,2) = 3
673  pos_trick(2,1) = 3
674  pos_trick(2,2) = 1
675  pos_trick(3,1) = 1
676  pos_trick(3,2) = 2
677  i2=pos_trick(i1, 1)
678  i3=pos_trick(i1, 2)
679  ip1=ine(i1, ie)
680  ip2=ine(i2, ie)
681  ip3=ine(i3, ie)
682  h=2.0*pdlib_tria(ie)
683  ugrad=-(y(ip3) - y(ip2))/h
684  vgrad= (x(ip3) - x(ip2))/h

References yowelementpool::ine, yownodepool::pdlib_tria, w3servmd::strace(), yownodepool::x, and yownodepool::y.

Referenced by trig_wave_setup_compute_system().

◆ trig_compute_lh_stress()

subroutine w3wavset::trig_compute_lh_stress ( real(rkind), dimension(npa), intent(out)  F_X,
real(rkind), dimension(npa), intent(out)  F_Y,
real(rkind), dimension(npa), intent(out)  DWNX 
)

Setup boundary pointer.

Parameters
[out]F_X
[out]F_Y
[out]DWNX
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
1-May-2018

Definition at line 449 of file w3wavset.F90.

449 !/
450 !/ +-----------------------------------+
451 !/ | WAVEWATCH III NOAA/NCEP |
452 !/ | |
453 !/ | Aron Roland (BGS IT&E GmbH) |
454 !/ | Mathieu Dutour-Sikiric (IRB) |
455 !/ | |
456 !/ | FORTRAN 90 |
457 !/ | Last update : 01-Mai-2018 |
458 !/ +-----------------------------------+
459 !/
460 !/ 01-Mai-2018 : Origination. ( version 6.04 )
461 !/
462 ! 1. Purpose : Setup boundary pointer
463 ! 2. Method :
464 ! 3. Parameters :
465 !
466 ! Parameter list
467 ! ----------------------------------------------------------------
468 ! ----------------------------------------------------------------
469 !
470 ! 4. Subroutines used :
471 !
472 ! Name Type Module Description
473 ! ----------------------------------------------------------------
474 ! STRACE Subr. W3SERVMD Subroutine tracing.
475 ! ----------------------------------------------------------------
476 !
477 ! 5. Called by :
478 !
479 ! Name Type Module Description
480 ! ----------------------------------------------------------------
481 ! ----------------------------------------------------------------
482 !
483 ! 6. Error messages :
484 ! 7. Remarks
485 ! 8. Structure :
486 ! 9. Switches :
487 !
488 ! !/S Enable subroutine tracing.
489 !
490 ! 10. Source code :
491 !
492 !/ ------------------------------------------------------------------- /
493 #ifdef W3_S
494  USE w3servmd, ONLY: strace
495 #endif
496 !
497  USE constants, ONLY: grav, dwat
498  use yownodepool, only: npa, iplg
499  USE w3gdatmd, only : mapfs
500  USE w3adatmd, ONLY: sxx, sxy, syy, wn, cg
501  USE w3parall, only: init_get_isea
502  USE w3odatmd, only : iaproc
503  USE w3gdatmd, ONLY : nseal, mapsta
504  USE w3adatmd, ONLY: dw
505  IMPLICIT NONE
506 !/
507 !/ ------------------------------------------------------------------- /
508 !/ Parameter list
509 !/
510 !/ ------------------------------------------------------------------- /
511 !/ Local PARAMETERs
512 !/
513 #ifdef W3_S
514  INTEGER, SAVE :: IENT = 0
515 #endif
516 !/
517 !/ ------------------------------------------------------------------- /
518 !/
519  real(rkind), intent(out) :: F_X(npa), F_Y(npa), DWNX(npa)
520  REAL(rkind) :: h
521  REAL(rkind) :: SXX_X, SXX_Y
522  REAL(rkind) :: SXY_X, SXY_Y
523  REAL(rkind) :: SYY_X, SYY_Y
524  INTEGER I, IP, IX
525  INTEGER JSEA, ISEA
526  real(rkind) :: U_X1(npa), U_Y1(npa)
527  real(rkind) :: U_X2(npa), U_Y2(npa)
528  real(rkind) :: SXX_p(npa), SXY_p(npa), SYY_p(npa)
529  real(rkind) :: eSXX, eSXY, eSYY
530  integer :: SXXmethod = 1
531 #ifdef W3_S
532  CALL strace (ient, 'VA_SETUP_IOBPD')
533 #endif
534  sxx_p=0
535  sxy_p=0
536  syy_p=0
537  dwnx=0
538  DO jsea=1,nseal
539  ip = jsea ! We remove the Z_status because now NX = NSEA
540  ix=iplg(ip)
541  isea=mapfs(1,ix)
542  IF (sxxmethod .eq. 1) THEN
543  esxx=sxx(jsea)/(dwat*grav)
544  esxy=sxy(jsea)/(dwat*grav)
545  esyy=syy(jsea)/(dwat*grav)
546  END IF
547  sxx_p(ip)=dble(esxx)
548  sxy_p(ip)=dble(esxy)
549  syy_p(ip)=dble(esyy)
550  dwnx(ip)=dw(isea)
551  END DO
552  !
553 #ifdef W3_DEBUGSTP
554  WRITE(740+iaproc,*) 'min/max(DEP)=', minval(dwnx), maxval(dwnx)
555  WRITE(740+iaproc,*) 'sum(abs(SXX))=', sum(abs(sxx_p))
556  WRITE(740+iaproc,*) 'sum(abs(SXY))=', sum(abs(sxy_p))
557  WRITE(740+iaproc,*) 'sum(abs(SYY))=', sum(abs(syy_p))
558  FLUSH(740+iaproc)
559 #endif
560 
561  CALL differentiate_xydir(sxx_p, u_x1, u_y1)
562 #ifdef W3_DEBUGSTP
563  WRITE(740+iaproc,*) 'sum(absU_XY1)=', sum(abs(u_x1)), sum(abs(u_y1))
564  FLUSH(740+iaproc)
565 #endif
566  CALL differentiate_xydir(sxy_p, u_x2, u_y2)
567 #ifdef W3_DEBUGSTP
568  WRITE(740+iaproc,*) 'sum(absU_XY2)=', sum(abs(u_x2)), sum(abs(u_y2))
569  FLUSH(740+iaproc)
570 #endif
571  f_x = -u_x1 - u_y2
572  !
573  CALL differentiate_xydir(syy_p, u_x1, u_y1)
574 #ifdef W3_DEBUGSTP
575  WRITE(740+iaproc,*) 'sum(absU_XY1)=', sum(abs(u_x1)), sum(abs(u_y1))
576  FLUSH(740+iaproc)
577 #endif
578  f_y = -u_y1 - u_x2
579 #ifdef W3_DEBUGSTP
580  WRITE(740+iaproc,*) 'sum(F_X)=', sum(f_x)
581  WRITE(740+iaproc,*) 'sum(F_Y)=', sum(f_y)
582  FLUSH(740+iaproc)
583 #endif

References w3adatmd::cg, differentiate_xydir(), w3adatmd::dw, constants::dwat, constants::grav, w3odatmd::iaproc, w3parall::init_get_isea(), yownodepool::iplg, w3gdatmd::mapfs, w3gdatmd::mapsta, yownodepool::npa, w3gdatmd::nseal, w3servmd::strace(), w3adatmd::sxx, w3adatmd::sxy, w3adatmd::syy, and w3adatmd::wn.

Referenced by trig_wave_setup_computation().

◆ trig_set_meanvalue_to_zero()

subroutine w3wavset::trig_set_meanvalue_to_zero ( real(rkind), dimension(npa), intent(inout)  TheVar)

Set mean value.

Parameters
[in,out]TheVar
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1380 of file w3wavset.F90.

1380 !/
1381 !/ +-----------------------------------+
1382 !/ | WAVEWATCH III NOAA/NCEP |
1383 !/ | |
1384 !/ | Mathieu Dutour-Sikiric (IRB) |
1385 !/ | Aron Roland (BGS IT&E GmbH) |
1386 !/ | |
1387 !/ | FORTRAN 90 |
1388 !/ | Last update : 01-Mai-2018 |
1389 !/ +-----------------------------------+
1390 !/
1391 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1392 !/
1393 ! 1. Purpose : set. mean value
1394 ! 2. Method :
1395 ! 3. Parameters :
1396 !
1397 ! Parameter list
1398 ! ----------------------------------------------------------------
1399 ! ----------------------------------------------------------------
1400 !
1401 ! 4. Subroutines used :
1402 !
1403 ! Name Type Module Description
1404 ! ----------------------------------------------------------------
1405 ! STRACE Subr. W3SERVMD Subroutine tracing.
1406 ! ----------------------------------------------------------------
1407 !
1408 ! 5. Called by :
1409 !
1410 ! Name Type Module Description
1411 ! ----------------------------------------------------------------
1412 ! ----------------------------------------------------------------
1413 !
1414 ! 6. Error messages :
1415 ! 7. Remarks
1416 ! 8. Structure :
1417 ! 9. Switches :
1418 !
1419 ! !/S Enable subroutine tracing.
1420 !
1421 ! 10. Source code :
1422 !
1423 !/ ------------------------------------------------------------------- /
1424 #ifdef W3_S
1425  USE w3servmd, ONLY: strace
1426 #endif
1427 !
1428  USE yownodepool, only: pdlib_si
1429  USE w3gdatmd, ONLY: nx, si
1430  USE w3gdatmd, ONLY: nseal
1431  USE w3adatmd, ONLY: mpi_comm_wcmp
1432  USE w3odatmd, only : iaproc, naproc, ntproc
1433  use yowdatapool, only: rtype, istatus
1434  use yownodepool, only: np, npa
1435  USE mpi, only : mpi_sum
1436  IMPLICIT NONE
1437 !/
1438 !/ ------------------------------------------------------------------- /
1439 !/ Parameter list
1440 !/
1441 !/ ------------------------------------------------------------------- /
1442 !/ Local PARAMETERs
1443 !/
1444 #ifdef W3_S
1445  INTEGER, SAVE :: IENT = 0
1446 #endif
1447 !/
1448 !/ ------------------------------------------------------------------- /
1449 !/
1450  real(rkind), intent(inout) :: TheVar(npa)
1451  real(rkind) :: SUM_SI_Var, SUM_SI, TheMean
1452  INTEGER IP, ierr
1453  real(rkind) :: eVect_loc(2), eVect_gl(2)
1454  integer iProc
1455 #ifdef W3_S
1456  CALL strace (ient, 'VA_SETUP_IOBPD')
1457 #endif
1458  sum_si_var=0
1459  sum_si=0
1460  DO ip=1,np
1461  sum_si_var = sum_si_var + pdlib_si(ip)*thevar(ip)
1462  sum_si = sum_si + pdlib_si(ip)
1463  END DO
1464  evect_loc(1)=sum_si_var
1465  evect_loc(2)=sum_si
1466 #ifdef W3_DEBUGSTP
1467  WRITE(740+iaproc,*) 'SUM_SI_Var=', sum_si_var, 'SUM_SI=', sum_si
1468  FLUSH(740+iaproc)
1469 #endif
1470  CALL mpi_allreduce(evect_loc,evect_gl,2,rtype,mpi_sum,mpi_comm_wcmp,ierr)
1471  sum_si_var=evect_gl(1)
1472  sum_si =evect_gl(2)
1473  themean=sum_si_var/sum_si
1474 #ifdef W3_DEBUGSTP
1475  WRITE(740+iaproc,*) 'TheMean=', themean
1476  FLUSH(740+iaproc)
1477 #endif
1478  DO ip=1,npa
1479  thevar(ip)=thevar(ip) - themean
1480  END DO

References w3odatmd::iaproc, yowdatapool::istatus, w3adatmd::mpi_comm_wcmp, w3odatmd::naproc, yownodepool::np, yownodepool::npa, w3gdatmd::nseal, w3odatmd::ntproc, w3gdatmd::nx, yownodepool::pdlib_si, yowdatapool::rtype, w3gdatmd::si, and w3servmd::strace().

Referenced by trig_wave_setup_computation().

◆ trig_wave_setup_apply_fct()

subroutine w3wavset::trig_wave_setup_apply_fct ( real(rkind), dimension(pdlib_nnz), intent(in)  ASPAR,
real(rkind), dimension(npa), intent(in)  TheIn,
real(rkind), dimension(npa), intent(out)  TheOut,
integer, dimension(npa), intent(in)  ACTIVE,
integer, dimension(npa), intent(in)  ACTIVESEC 
)
Parameters
[in]ASPAR
[in]TheIn
[out]TheOut
[in]ACTIVE
[in]ACTIVESEC
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1021 of file w3wavset.F90.

1021 !/
1022 !/ +-----------------------------------+
1023 !/ | WAVEWATCH III NOAA/NCEP |
1024 !/ | |
1025 !/ | Mathieu Dutour-Sikiric (IRB) |
1026 !/ | Aron Roland (BGS IT&E GmbH) |
1027 !/ | |
1028 !/ | FORTRAN 90 |
1029 !/ | Last update : 01-Mai-2018 |
1030 !/ +-----------------------------------+
1031 !/
1032 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1033 !/
1034 ! 1. Purpose : compute off diagonal contr.
1035 ! 2. Method :
1036 ! 3. Parameters :
1037 !
1038 ! Parameter list
1039 ! ----------------------------------------------------------------
1040 ! ----------------------------------------------------------------
1041 !
1042 ! 4. Subroutines used :
1043 !
1044 ! Name Type Module Description
1045 ! ----------------------------------------------------------------
1046 ! STRACE Subr. W3SERVMD Subroutine tracing.
1047 ! ----------------------------------------------------------------
1048 !
1049 ! 5. Called by :
1050 !
1051 ! Name Type Module Description
1052 ! ----------------------------------------------------------------
1053 ! ----------------------------------------------------------------
1054 !
1055 ! 6. Error messages :
1056 ! 7. Remarks
1057 ! 8. Structure :
1058 ! 9. Switches :
1059 !
1060 ! !/S Enable subroutine tracing.
1061 !
1062 ! 10. Source code :
1063 !
1064 !/ ------------------------------------------------------------------- /
1065 #ifdef W3_S
1066  USE w3servmd, ONLY: strace
1067 #endif
1068 !
1070  USE yownodepool, only: pdlib_ia, pdlib_ja, pdlib_nnz
1071  use yownodepool, only: np, npa
1072  USE w3gdatmd, ONLY: nseal
1073  IMPLICIT NONE
1074 !/
1075 !/ ------------------------------------------------------------------- /
1076 !/ Parameter list
1077 !/
1078 !/ ------------------------------------------------------------------- /
1079 !/ Local PARAMETERs
1080 !/
1081 #ifdef W3_S
1082  INTEGER, SAVE :: IENT = 0
1083 #endif
1084 !/
1085 !/ ------------------------------------------------------------------- /
1086 !/
1087  REAL(rkind), intent(in) :: ASPAR(PDLIB_NNZ)
1088  REAL(rkind), intent(in) :: TheIn(npa)
1089  REAL(rkind), intent(out) :: TheOut(npa)
1090  INTEGER, intent(in) :: ACTIVE(npa), ACTIVESEC(npa)
1091  integer IP, J, JP
1092  REAL(rkind) :: eCoeff
1093 #ifdef W3_S
1094  CALL strace (ient, 'VA_SETUP_IOBPD')
1095 #endif
1096  theout=0
1097  DO ip=1,npa
1098  IF (activesec(ip) .eq. 1) THEN
1099  DO j=pdlib_ia(ip),pdlib_ia(ip+1)-1
1100  jp=pdlib_ja(j)
1101  ecoeff=aspar(j)
1102  theout(ip)=theout(ip) + ecoeff*thein(jp)
1103  END DO
1104  END IF
1105  END DO
1106  CALL pdlib_exchange1dreal(theout)

References yownodepool::np, yownodepool::npa, w3gdatmd::nseal, yowexchangemodule::pdlib_exchange1dreal(), yownodepool::pdlib_ia, yownodepool::pdlib_ja, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by trig_wave_setup_solve_poisson_neumann_dir().

◆ trig_wave_setup_apply_precond()

subroutine w3wavset::trig_wave_setup_apply_precond ( real(rkind), dimension(pdlib_nnz), intent(in)  ASPAR,
real(rkind), dimension(npa), intent(in)  TheIn,
real(rkind), dimension(npa), intent(out)  TheOut,
integer, dimension(npa), intent(in)  ACTIVE,
integer, dimension(npa), intent(in)  ACTIVESEC 
)

Preconditioner.

Parameters
[in]ASPAR
[in]TheIn
[out]TheOut
[in]ACTIVE
[in]ACTIVESEC
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 893 of file w3wavset.F90.

893 !/
894 !/ +-----------------------------------+
895 !/ | WAVEWATCH III NOAA/NCEP |
896 !/ | |
897 !/ | Mathieu Dutour-Sikiric (IRB) |
898 !/ | Aron Roland (BGS IT&E GmbH) |
899 !/ | |
900 !/ | FORTRAN 90 |
901 !/ | Last update : 01-Mai-2018 |
902 !/ +-----------------------------------+
903 !/
904 !/ 01-Mai-2018 : Origination. ( version 6.04 )
905 !/
906 ! 1. Purpose : preconditioner
907 ! 2. Method :
908 ! 3. Parameters :
909 !
910 ! Parameter list
911 ! ----------------------------------------------------------------
912 ! ----------------------------------------------------------------
913 !
914 ! 4. Subroutines used :
915 !
916 ! Name Type Module Description
917 ! ----------------------------------------------------------------
918 ! STRACE Subr. W3SERVMD Subroutine tracing.
919 ! ----------------------------------------------------------------
920 !
921 ! 5. Called by :
922 !
923 ! Name Type Module Description
924 ! ----------------------------------------------------------------
925 ! ----------------------------------------------------------------
926 !
927 ! 6. Error messages :
928 ! 7. Remarks
929 ! 8. Structure :
930 ! 9. Switches :
931 !
932 ! !/S Enable subroutine tracing.
933 !
934 ! 10. Source code :
935 !
936 !/ ------------------------------------------------------------------- /
937 #ifdef W3_S
938  USE w3servmd, ONLY: strace
939 #endif
940 !
943  use yownodepool, only: npa
944  USE w3odatmd, only : iaproc
945  USE w3odatmd, only : iaproc
946  USE yownodepool, only: iplg
947  IMPLICIT NONE
948 !/
949 !/ ------------------------------------------------------------------- /
950 !/ Parameter list
951 !/
952 !/ ------------------------------------------------------------------- /
953 !/ Local PARAMETERs
954 !/
955 #ifdef W3_S
956  INTEGER, SAVE :: IENT = 0
957 #endif
958 !/
959 !/ ------------------------------------------------------------------- /
960 !/
961  REAL(rkind), intent(in) :: ASPAR(PDLIB_NNZ)
962  REAL(rkind), intent(in) :: TheIn(npa)
963  REAL(rkind), intent(out) :: TheOut(npa)
964  INTEGER, intent(IN) :: ACTIVE(npa), ACTIVESEC(npa)
965  integer IP, J1, J, JP, J2
966  REAL(rkind) :: eCoeff
967  INTEGER :: ThePrecond = 2
968 #ifdef W3_S
969  CALL strace (ient, 'VA_SETUP_IOBPD')
970 #endif
971  IF (theprecond .eq. 0) THEN
972  theout=thein
973  END IF
974  IF (theprecond .eq. 1) THEN
975  theout=0
976  DO ip=1,npa
977  IF (active(ip) .eq. 1) THEN
978  j1=pdlib_i_diag(ip)
979  DO j=pdlib_ia(ip),pdlib_ia(ip+1)-1
980  jp=pdlib_ja(j)
981  IF (activesec(jp) .eq. 1) THEN
982  IF (j .eq. j1) THEN
983  ecoeff=1.0/aspar(j)
984  ELSE
985  j2=pdlib_i_diag(jp)
986  ecoeff=-aspar(j) /(aspar(j1)*aspar(j2))
987  END IF
988  theout(ip)=theout(ip) + ecoeff*thein(jp)
989  END IF
990  END DO
991  END IF
992  END DO
993  END IF
994  IF (theprecond .eq. 2) THEN
995  DO ip=1,npa
996  IF (activesec(ip) .eq. 1) THEN
997  j=pdlib_i_diag(ip)
998  theout(ip)=thein(ip)/aspar(j)
999  ELSE
1000  theout(ip)=thein(ip)
1001  END IF
1002  END DO
1003  END IF
1004  CALL pdlib_exchange1dreal(theout)

References w3odatmd::iaproc, yownodepool::iplg, yownodepool::npa, yowexchangemodule::pdlib_exchange1dreal(), yownodepool::pdlib_i_diag, yownodepool::pdlib_ia, yownodepool::pdlib_ja, yownodepool::pdlib_nnz, and w3servmd::strace().

Referenced by trig_wave_setup_solve_poisson_neumann_dir().

◆ trig_wave_setup_computation()

subroutine w3wavset::trig_wave_setup_computation

Setup computation.

Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1598 of file w3wavset.F90.

1598 !/
1599 !/ +-----------------------------------+
1600 !/ | WAVEWATCH III NOAA/NCEP |
1601 !/ | |
1602 !/ | Mathieu Dutour-Sikiric (IRB) |
1603 !/ | Aron Roland (BGS IT&E GmbH) |
1604 !/ | |
1605 !/ | FORTRAN 90 |
1606 !/ | Last update : 01-Mai-2018 |
1607 !/ +-----------------------------------+
1608 !/
1609 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1610 !/
1611 ! 1. Purpose : Setup computation
1612 ! 2. Method :
1613 ! 3. Parameters :
1614 !
1615 ! Parameter list
1616 ! ----------------------------------------------------------------
1617 ! ----------------------------------------------------------------
1618 !
1619 ! 4. Subroutines used :
1620 !
1621 ! Name Type Module Description
1622 ! ----------------------------------------------------------------
1623 ! STRACE Subr. W3SERVMD Subroutine tracing.
1624 ! ----------------------------------------------------------------
1625 !
1626 ! 5. Called by :
1627 !
1628 ! Name Type Module Description
1629 ! ----------------------------------------------------------------
1630 ! ----------------------------------------------------------------
1631 !
1632 ! 6. Error messages :
1633 ! 7. Remarks
1634 ! 8. Structure :
1635 ! 9. Switches :
1636 !
1637 ! !/S Enable subroutine tracing.
1638 !
1639 ! 10. Source code :
1640 !
1641 !/ ------------------------------------------------------------------- /
1642 #ifdef W3_S
1643  USE w3servmd, ONLY: strace
1644 #endif
1645 !
1646  USE yownodepool, only: pdlib_nnz, pdlib_ia, pdlib_ja, iplg, npa, np
1647  USE w3gdatmd, only : mapfs
1649  USE w3adatmd, ONLY: dw
1650  USE w3gdatmd, ONLY: nseal, nsea, nx
1651  USE w3wdatmd, ONLY: zeta_setup
1652  USE w3odatmd, only : iaproc, naproc, ntproc
1653  USE w3parall, only: init_get_isea
1655  IMPLICIT NONE
1656 !/
1657 !/ ------------------------------------------------------------------- /
1658 !/ Parameter list
1659 !/
1660 !/ ------------------------------------------------------------------- /
1661 !/ Local PARAMETERs
1662 !/
1663 #ifdef W3_S
1664  INTEGER, SAVE :: IENT = 0
1665 #endif
1666 !/
1667 !/ ------------------------------------------------------------------- /
1668 !/
1669 !
1670 ! CALL W3SETG
1671  REAL(rkind) :: ZETA_WORK(npa)
1672  REAL(rkind) :: ZETA_WORK_ALL(NX)
1673  REAL(rkind) :: F_X(npa), F_Y(npa), DWNX(npa)
1674  REAL(rkind) :: ASPAR(PDLIB_NNZ), B(npa)
1675  INTEGER I, ISEA, JSEA, IX, IP, IP_glob
1676  INTEGER :: ACTIVE(npa), ACTIVESEC(npa)
1677  REAL(rkind) max_val, min_val
1678 #ifdef W3_S
1679  CALL strace (ient, 'VA_SETUP_IOBPD')
1680 #endif
1681 ! ZETA_SETUP is allocated on 1:NSEA
1682 ! ZETA_WORK is on 1:npa
1683 #ifdef W3_DEBUGSTP
1684  WRITE(740+iaproc,*) 'NAPROC=', naproc, ' NTPROC=', ntproc
1685  WRITE(740+iaproc,*) 'NSEAL=', nseal
1686  WRITE(740+iaproc,*) 'npa=', npa, ' np=', np
1687  FLUSH(740+iaproc)
1688 #endif
1689  zeta_work=0
1690  DO ip=1,npa
1691  ix=iplg(ip)
1692  isea=mapfs(1,ix)
1693  IF (isea .gt. 0) THEN
1694  zeta_work(ip)=zeta_setup(isea)
1695  END IF
1696  END DO
1697 #ifdef W3_DEBUGSTP
1698  WRITE(740+iaproc,*) 'Before TRIG_COMPUTE_LH_STRESS'
1699  FLUSH(740+iaproc)
1700 #endif
1701 
1702  CALL trig_compute_lh_stress(f_x, f_y, dwnx)
1703 #ifdef W3_DEBUGSTP
1704  WRITE(740+iaproc,*) 'After TRIG_COMPUTE_LH_STRESS'
1705  FLUSH(740+iaproc)
1706 #endif
1707  CALL compute_active_node(dwnx, active)
1708 #ifdef W3_DEBUGSTP
1709  WRITE(740+iaproc,*) 'After COMPUTE_ACTIVE_NODE'
1710  FLUSH(740+iaproc)
1711 #endif
1712  CALL trig_wave_setup_compute_system(aspar, b, f_x, f_y, dwnx, active, activesec)
1713 #ifdef W3_DEBUGSTP
1714  WRITE(740+iaproc,*) 'Before,B,min=', minval(b), ' max=', maxval(b)
1715  FLUSH(740+iaproc)
1716 #endif
1717 
1718 
1719 ! CALL TRIG_SET_MEANVALUE_TO_ZERO(B)
1720 #ifdef W3_DEBUGSTP
1721  WRITE(740+iaproc,*) 'After,B,min=', minval(b), ' max=', maxval(b)
1722  FLUSH(740+iaproc)
1723 #endif
1724 
1725 
1726  CALL trig_wave_setup_solve_poisson_neumann_dir(aspar, b, zeta_work, active, activesec)
1727 
1728  CALL trig_set_meanvalue_to_zero(zeta_work)
1729 #ifdef W3_DEBUGSTP
1730  WRITE(740+iaproc,*) 'After SET_MEAN ZETA_WORK(min/max)=', minval(zeta_work), maxval(zeta_work)
1731  FLUSH(740+iaproc)
1732 #endif
1733  CALL pdlib_exchange1dreal(zeta_work)
1734  max_val = -100000000
1735  min_val = -100000000
1736  DO ip=1,npa
1737  ix=iplg(ip)
1738  isea=mapfs(1,ix)
1739  IF (isea .gt. 0) THEN
1740  zeta_setup(isea) = zeta_work(ip)
1741  max_val = max(max_val, zeta_work(ip))
1742  min_val = max(min_val, zeta_work(ip))
1743  END IF
1744  END DO
1745 #ifdef W3_DEBUGSTP
1746  WRITE(740+iaproc,*) 'TRIG_WAVE_SETUP_COMPUTATION, max/min=', max_val, min_val
1747  FLUSH(740+iaproc)
1748 #endif
1749  zeta_work_all = 0.
1750  DO ip = 1, npa
1751  isea = iplg(ip)
1752  zeta_work_all(isea) = zeta_work(ip)
1753  END DO
1754  CALL synchronize_global_array(zeta_work_all)
1755  DO ix = 1, nx
1756  zeta_setup(ix) = zeta_work_all(ix)
1757  END DO
1758  IF (iaproc .EQ. 1) THEN
1759  write(6666) 1.
1760  write(6666) (zeta_work_all(ix), zeta_work_all(ix), zeta_work_all(ix), ix = 1, nx)
1761  ENDIF
1762 #ifdef W3_DEBUGSTP
1763  WRITE(740+iaproc,*) 'Now exiting TRIG_WAVE_SETUP_COMPUTATION'
1764  FLUSH(740+iaproc)
1765 #endif

References compute_active_node(), w3adatmd::dw, w3odatmd::iaproc, w3parall::init_get_isea(), yownodepool::iplg, w3gdatmd::mapfs, w3odatmd::naproc, yownodepool::np, yownodepool::npa, w3gdatmd::nsea, w3gdatmd::nseal, w3odatmd::ntproc, w3gdatmd::nx, yowexchangemodule::pdlib_exchange1dreal(), yownodepool::pdlib_ia, yownodepool::pdlib_ja, yownodepool::pdlib_nnz, w3servmd::strace(), w3parall::synchronize_global_array(), trig_compute_lh_stress(), trig_set_meanvalue_to_zero(), trig_wave_setup_compute_system(), trig_wave_setup_solve_poisson_neumann_dir(), and w3wdatmd::zeta_setup.

Referenced by wave_setup_computation().

◆ trig_wave_setup_compute_system()

subroutine w3wavset::trig_wave_setup_compute_system ( real(rkind), dimension(pdlib_nnz), intent(out)  ASPAR,
real(rkind), dimension(npa), intent(out)  B,
real(rkind), dimension(npa), intent(in)  FX,
real(rkind), dimension(npa), intent(in)  FY,
real(rkind), dimension(npa), intent(in)  DWNX,
integer, dimension(npa), intent(in)  ACTIVE,
integer, dimension(npa), intent(out)  ACTIVESEC 
)

Setup system matrix for solutions of wave setup eq.

Parameters
[in]FX
[in]FY
[in]DWNX
[out]ASPAR
[out]B
[in]ACTIVE
[out]ACTIVESEC
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 703 of file w3wavset.F90.

703 !/
704 !/ +-----------------------------------+
705 !/ | WAVEWATCH III NOAA/NCEP |
706 !/ | |
707 !/ | Mathieu Dutour-Sikiric (IRB) |
708 !/ | Aron Roland (BGS IT&E GmbH) |
709 !/ | |
710 !/ | FORTRAN 90 |
711 !/ | Last update : 01-Mai-2018 |
712 !/ +-----------------------------------+
713 !/
714 !/ 01-Mai-2018 : Origination. ( version 6.04 )
715 !/
716 ! 1. Purpose : Setup system matrix for solutions of wave setup eq.
717 ! 2. Method :
718 ! 3. Parameters :
719 !
720 ! Parameter list
721 ! ----------------------------------------------------------------
722 ! ----------------------------------------------------------------
723 !
724 ! 4. Subroutines used :
725 !
726 ! Name Type Module Description
727 ! ----------------------------------------------------------------
728 ! STRACE Subr. W3SERVMD Subroutine tracing.
729 ! ----------------------------------------------------------------
730 !
731 ! 5. Called by :
732 !
733 ! Name Type Module Description
734 ! ----------------------------------------------------------------
735 ! ----------------------------------------------------------------
736 !
737 ! 6. Error messages :
738 ! 7. Remarks
739 ! 8. Structure :
740 ! 9. Switches :
741 !
742 ! !/S Enable subroutine tracing.
743 !
744 ! 10. Source code :
745 !
746 !/ ------------------------------------------------------------------- /
747 #ifdef W3_S
748  USE w3servmd, ONLY: strace
749 #endif
750 !
751  use yowelementpool, only: ine, ne
753  use yownodepool, only: pdlib_i_diag
754  USE yownodepool, only: iplg
755  USE w3odatmd, only : iaproc
756  IMPLICIT NONE
757 !/
758 !/ ------------------------------------------------------------------- /
759 !/ Parameter list
760 !/
761 !/ ------------------------------------------------------------------- /
762 !/ Local PARAMETERs
763 !/
764 #ifdef W3_S
765  INTEGER, SAVE :: IENT = 0
766 #endif
767 !/
768 !/ ------------------------------------------------------------------- /
769 !/
770  real(rkind), intent(in) :: FX(npa), FY(npa), DWNX(npa)
771  real(rkind), intent(out) :: ASPAR(PDLIB_NNZ)
772  real(rkind), intent(out) :: B(npa)
773  integer, intent(in) :: ACTIVE(npa)
774  integer, intent(out) :: ACTIVESEC(npa)
775  INTEGER :: POS_TRICK(3,2), POS_SHIFT(3,3)
776  integer I1, I2, I3, IP1, IP2, IP3
777  integer IDX, IDX1, IDX2, IDX3
778  INTEGER IE, IP, I, J, K, IPp, JPp
779  real(rkind) :: eDep, eFX, eFY, eScal, eFact, eArea
780  real(rkind) :: UGRAD, VGRAD, UGRAD1, VGRAD1
781  real(rkind) :: eOff
782  logical DoPrintOut
783  INTEGER sumActive
784  INTEGER LIDX(2), KIDX(2), jdx
785  INTEGER IPglob1, IPglob2, IPglob3
786 #ifdef W3_DEBUGSTP
787  REAL(rkind) :: ListDiag(npa)
788 #endif
789 #ifdef W3_S
790  CALL strace (ient, 'VA_SETUP_IOBPD')
791 #endif
792  pos_trick(1,1) = 2
793  pos_trick(1,2) = 3
794  pos_trick(2,1) = 3
795  pos_trick(2,2) = 1
796  pos_trick(3,1) = 1
797  pos_trick(3,2) = 2
798  aspar=0
799  b=0
800  DO i=1,3
801  DO j=1,3
802  k= i-j+1
803  IF (k .le. 0) THEN
804  k=k+3
805  END IF
806  IF (k .ge. 4) THEN
807  k=k-3
808  END IF
809  pos_shift(i,j)=k
810  END DO
811  END DO
812  DO i=1,3
813  jdx=0
814  DO idx=1,3
815  k=pos_shift(i,idx)
816  IF (k .ne. i) THEN
817  jdx=jdx+1
818  lidx(jdx)=idx
819  kidx(jdx)=k
820  END IF
821  END DO
822  pos_shift(i,lidx(1))=kidx(2)
823  pos_shift(i,lidx(2))=kidx(1)
824  END DO
825  activesec=0
826  DO ie=1,ne
827  ip1=ine(1,ie)
828  ip2=ine(2,ie)
829  ip3=ine(3,ie)
830  efx =(fx(ip1) + fx(ip2) + fx(ip3))/3
831  efy =(fy(ip1) + fy(ip2) + fy(ip3))/3
832  sumactive=active(ip1) + active(ip2) + active(ip3)
833  IF (sumactive .eq. 3) THEN
834  activesec(ip1)=1
835  activesec(ip2)=1
836  activesec(ip3)=1
837  edep=(dwnx(ip1) + dwnx(ip2) + dwnx(ip3))/3.0
838  earea=pdlib_tria(ie)
839  efact=edep*earea
840  DO i1=1,3
841  i2=pos_trick(i1,1)
842  i3=pos_trick(i1,2)
843  ip1=ine(i1,ie)
844  ip2=ine(i2,ie)
845  ip3=ine(i3,ie)
846  idx1=pdlib_ja_ie(i1,1,ie)
847  idx2=pdlib_ja_ie(i1,2,ie)
848  idx3=pdlib_ja_ie(i1,3,ie)
849  CALL trig_compute_diff(ie, i1, ugrad1, vgrad1)
850  escal=ugrad1*efx + vgrad1*efy
851  b(ip1) = b(ip1) + escal*earea
852  !
853  DO idx=1,3
854  k=pos_shift(i1, idx)
855  CALL trig_compute_diff(ie, k, ugrad, vgrad)
856  escal=ugrad*ugrad1 + vgrad*vgrad1
857  j=pdlib_ja_ie(i1,idx,ie)
858  aspar(j)=aspar(j) + efact*escal
859  END DO
860  END DO
861  END IF
862  END DO
863  doprintout=.true.
864  IF (doprintout .eqv. .true.) THEN
865  DO ip=1,np
866  eoff=0
867  END DO
868  END IF
869 #ifdef W3_DEBUGSTP
870  DO ip=1,npa
871  j=pdlib_i_diag(ip)
872  listdiag(ip)=aspar(j)
873  END DO
874  WRITE(740+iaproc,*) 'Diag, min=', minval(listdiag), ' max=', maxval(listdiag)
875  WRITE(740+iaproc,*) 'Diag, quot=', maxval(listdiag)/minval(listdiag)
876 #endif

References w3odatmd::iaproc, yowelementpool::ine, yownodepool::iplg, yowelementpool::ne, yownodepool::np, yownodepool::npa, yownodepool::pdlib_i_diag, yownodepool::pdlib_ja_ie, yownodepool::pdlib_nnz, yownodepool::pdlib_tria, w3servmd::strace(), and trig_compute_diff().

Referenced by trig_wave_setup_computation().

◆ trig_wave_setup_scalar_prod()

subroutine w3wavset::trig_wave_setup_scalar_prod ( real(rkind), dimension(npa), intent(in)  V1,
real(rkind), dimension(npa), intent(in)  V2,
real(rkind), intent(inout)  eScal 
)

Scalar product plus exchange.

Parameters
[in]V1
[in]V2
[in,out]eScal
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1121 of file w3wavset.F90.

1121 !/
1122 !/ +-----------------------------------+
1123 !/ | WAVEWATCH III NOAA/NCEP |
1124 !/ | |
1125 !/ | Mathieu Dutour-Sikiric (IRB) |
1126 !/ | Aron Roland (BGS IT&E GmbH) |
1127 !/ | |
1128 !/ | FORTRAN 90 |
1129 !/ | Last update : 01-Mai-2018 |
1130 !/ +-----------------------------------+
1131 !/
1132 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1133 !/
1134 ! 1. Purpose : scalar prod. + exchange
1135 ! 2. Method :
1136 ! 3. Parameters :
1137 !
1138 ! Parameter list
1139 ! ----------------------------------------------------------------
1140 ! ----------------------------------------------------------------
1141 !
1142 ! 4. Subroutines used :
1143 !
1144 ! Name Type Module Description
1145 ! ----------------------------------------------------------------
1146 ! STRACE Subr. W3SERVMD Subroutine tracing.
1147 ! ----------------------------------------------------------------
1148 !
1149 ! 5. Called by :
1150 !
1151 ! Name Type Module Description
1152 ! ----------------------------------------------------------------
1153 ! ----------------------------------------------------------------
1154 !
1155 ! 6. Error messages :
1156 ! 7. Remarks
1157 ! 8. Structure :
1158 ! 9. Switches :
1159 !
1160 ! !/S Enable subroutine tracing.
1161 !
1162 ! 10. Source code :
1163 !
1164 !/ ------------------------------------------------------------------- /
1165 #ifdef W3_S
1166  USE w3servmd, ONLY: strace
1167 #endif
1168 !
1169  USE w3gdatmd, ONLY: nx
1170  USE w3adatmd, ONLY: mpi_comm_wcmp
1171  use yowdatapool, only: rtype, istatus
1172  use yownodepool, only: np, npa
1173  USE w3odatmd, only : iaproc, naproc, ntproc
1174  USE w3gdatmd, ONLY: nseal
1175  USE mpi, only : mpi_sum
1176  IMPLICIT NONE
1177 !/
1178 !/ ------------------------------------------------------------------- /
1179 !/ Parameter list
1180 !/
1181 !/ ------------------------------------------------------------------- /
1182 !/ Local PARAMETERs
1183 !/
1184 #ifdef W3_S
1185  INTEGER, SAVE :: IENT = 0
1186 #endif
1187 !/
1188 !/ ------------------------------------------------------------------- /
1189 !/
1190  real(rkind), intent(in) :: V1(npa), V2(npa)
1191  real(rkind), intent(inout) :: eScal
1192  integer IP
1193  real(rkind) :: lScal_loc(1), lScal_gl(1)
1194  integer ierr
1195 #ifdef W3_S
1196  CALL strace (ient, 'VA_SETUP_IOBPD')
1197 #endif
1198  lscal_loc = 0
1199  DO ip=1,np
1200  lscal_loc(1) = lscal_loc(1) + v1(ip)*v2(ip)
1201  END DO
1202  CALL mpi_allreduce(lscal_loc,lscal_gl,1,rtype,mpi_sum,mpi_comm_wcmp,ierr)
1203  escal = lscal_gl(1)

References w3odatmd::iaproc, yowdatapool::istatus, w3adatmd::mpi_comm_wcmp, w3odatmd::naproc, yownodepool::np, yownodepool::npa, w3gdatmd::nseal, w3odatmd::ntproc, w3gdatmd::nx, yowdatapool::rtype, and w3servmd::strace().

Referenced by trig_wave_setup_solve_poisson_neumann_dir().

◆ trig_wave_setup_solve_poisson_neumann_dir()

subroutine w3wavset::trig_wave_setup_solve_poisson_neumann_dir ( real(rkind), dimension(pdlib_nnz), intent(in)  ASPAR,
real(rkind), dimension(npa), intent(in)  B,
real(rkind), dimension(npa), intent(out)  TheOut,
integer, dimension(npa), intent(in)  ACTIVE,
integer, dimension(npa), intent(in)  ACTIVESEC 
)

Poisson equation solver.

Parameters
[in]ASPAR
[in]B
[out]TheOut
[in]ACTIVE
[in]ACTIVESEC
Author
Mathieu Dutour-Sikiric
Aron Roland
Date
1-May-2018

Definition at line 1220 of file w3wavset.F90.

1220 !/
1221 !/ +-----------------------------------+
1222 !/ | WAVEWATCH III NOAA/NCEP |
1223 !/ | |
1224 !/ | Mathieu Dutour-Sikiric (IRB) |
1225 !/ | Aron Roland (BGS IT&E GmbH) |
1226 !/ | |
1227 !/ | FORTRAN 90 |
1228 !/ | Last update : 01-Mai-2018 |
1229 !/ +-----------------------------------+
1230 !/
1231 !/ 01-Mai-2018 : Origination. ( version 6.04 )
1232 !/
1233 ! 1. Purpose : poisson eq. solver
1234 ! 2. Method :
1235 ! 3. Parameters :
1236 !
1237 ! Parameter list
1238 ! ----------------------------------------------------------------
1239 ! ----------------------------------------------------------------
1240 !
1241 ! 4. Subroutines used :
1242 !
1243 ! Name Type Module Description
1244 ! ----------------------------------------------------------------
1245 ! STRACE Subr. W3SERVMD Subroutine tracing.
1246 ! ----------------------------------------------------------------
1247 !
1248 ! 5. Called by :
1249 !
1250 ! Name Type Module Description
1251 ! ----------------------------------------------------------------
1252 ! ----------------------------------------------------------------
1253 !
1254 ! 6. Error messages :
1255 ! 7. Remarks
1256 ! 8. Structure :
1257 ! 9. Switches :
1258 !
1259 ! !/S Enable subroutine tracing.
1260 !
1261 ! 10. Source code :
1262 !
1263 !/ ------------------------------------------------------------------- /
1264 #ifdef W3_S
1265  USE w3servmd, ONLY: strace
1266 #endif
1267 !
1268  USE yownodepool, only: pdlib_nnz
1269  USE w3gdatmd, ONLY: nseal, solverthr_stp
1270  USE w3odatmd, only : iaproc
1271  use yownodepool, only: np, npa
1272  IMPLICIT NONE
1273 !/
1274 !/ ------------------------------------------------------------------- /
1275 !/ Parameter list
1276 !/
1277 !/ ------------------------------------------------------------------- /
1278 !/ Local PARAMETERs
1279 !/
1280 #ifdef W3_S
1281  INTEGER, SAVE :: IENT = 0
1282 #endif
1283 !/
1284 !/ ------------------------------------------------------------------- /
1285 !/
1286  real(rkind), intent(in) :: ASPAR(PDLIB_NNZ)
1287  real(rkind), intent(in) :: B(npa)
1288  real(rkind), intent(out) :: TheOut(npa)
1289  integer, intent(in) :: ACTIVE(npa), ACTIVESEC(npa)
1290  real(rkind) :: V_X(npa), V_R(npa), V_Z(npa), V_P(npa), V_Y(npa)
1291  real(rkind) :: uO, uN, alphaV, h1, h2
1292  real(rkind) :: eNorm, beta
1293  real(rkind) :: SOLVERTHR
1294  integer IP, nbIter
1295 #ifdef W3_S
1296  CALL strace (ient, 'VA_SETUP_IOBPD')
1297 #endif
1298  solverthr = solverthr_stp
1299 
1300 #ifdef W3_DEBUGSTP
1301  WRITE(740+iaproc,*) 'Begin TRIG_WAVE_SETUP_SOLVE ....'
1302  FLUSH(740+iaproc)
1303 #endif
1304  nbiter=0
1305  v_x=0
1306  v_r=b
1307  CALL trig_wave_setup_apply_precond(aspar, v_r, v_z, active, activesec)
1308  v_p=v_z
1309  CALL trig_wave_setup_scalar_prod(v_z, v_r, uo)
1310 #ifdef W3_DEBUGSTP
1311  WRITE(740+iaproc,*) 'uO=', uo
1312  FLUSH(740+iaproc)
1313 #endif
1314  CALL trig_wave_setup_scalar_prod(b, b, enorm)
1315 #ifdef W3_DEBUGSTP
1316  WRITE(740+iaproc,*) 'eNorm(B)=', enorm
1317  WRITE(740+iaproc,*) 'SOLVERTHR=', solverthr
1318  WRITE(740+iaproc,*) 'SOLVERTHR=', solverthr, ' eNorm(B)=', enorm
1319  FLUSH(740+iaproc)
1320 #endif
1321  IF (enorm .le. solverthr) THEN
1322 #ifdef W3_DEBUGSTP
1323  WRITE(740+iaproc,*) 'Leaving here, zero solution'
1324  FLUSH(740+iaproc)
1325 #endif
1326  theout=v_x
1327  RETURN
1328  END IF
1329  DO
1330  nbiter=nbiter + 1
1331  CALL trig_wave_setup_apply_fct(aspar, v_p, v_y, active, activesec)
1332  CALL trig_wave_setup_scalar_prod(v_p, v_y, h2)
1333  alphav=uo/h2
1334  !
1335  DO ip=1,npa
1336  v_x(ip) = v_x(ip) + alphav * v_p(ip)
1337  v_r(ip) = v_r(ip) - alphav * v_y(ip)
1338  END DO
1339  !
1340  CALL trig_wave_setup_scalar_prod(v_r, v_r, enorm)
1341 #ifdef W3_DEBUGSTP
1342  WRITE(740+iaproc,*) 'nbIter=', nbiter, ' eNorm(res)=', enorm
1343  FLUSH(740+iaproc)
1344 #endif
1345  IF (enorm .le. solverthr) THEN
1346  EXIT
1347  END IF
1348  !
1349  CALL trig_wave_setup_apply_precond(aspar, v_r, v_z, active, activesec)
1350  CALL trig_wave_setup_scalar_prod(v_z, v_r, un)
1351  !
1352  beta=un/uo
1353  uo=un
1354 #ifdef W3_DEBUGSTP
1355  WRITE(740+iaproc,*) ' beta=', beta, ' uN=', un, ' alphaV=', alphav, ' h2=', h2
1356  FLUSH(740+iaproc)
1357 #endif
1358  !
1359  DO ip=1,npa
1360  v_p(ip)=v_z(ip) + beta * v_p(ip)
1361  END DO
1362  END DO
1363  theout=v_x
1364 #ifdef W3_DEBUGSTP
1365  WRITE(740+iaproc,*) 'TRIG_WAVE_SETUP_SOLVE_POISSON_NEUMANN_DIR, max/min=', maxval(theout), minval(theout)
1366  FLUSH(740+iaproc)
1367 #endif

References w3odatmd::iaproc, yownodepool::np, yownodepool::npa, w3gdatmd::nseal, yownodepool::pdlib_nnz, w3gdatmd::solverthr_stp, w3servmd::strace(), trig_wave_setup_apply_fct(), trig_wave_setup_apply_precond(), and trig_wave_setup_scalar_prod().

Referenced by trig_wave_setup_computation().

◆ wave_setup_computation()

subroutine w3wavset::wave_setup_computation

General driver.

Author
Aron Roland
Mathieu Dutour-Sikiric
Date
1-May-2018

Definition at line 3075 of file w3wavset.F90.

3075 !/
3076 !/ +-----------------------------------+
3077 !/ | WAVEWATCH III NOAA/NCEP |
3078 !/ | |
3079 !/ | Aron Roland (BGS IT&E GmbH) |
3080 !/ | Mathieu Dutour-Sikiric (IRB) |
3081 !/ | |
3082 !/ | FORTRAN 90 |
3083 !/ | Last update : 01-Mai-2018 |
3084 !/ +-----------------------------------+
3085 !/
3086 !/ 01-Mai-2018 : Origination. ( version 6.04 )
3087 !/
3088 ! 1. Purpose : general driver
3089 ! 2. Method :
3090 ! 3. Parameters :
3091 !
3092 ! Parameter list
3093 ! ----------------------------------------------------------------
3094 ! ----------------------------------------------------------------
3095 !
3096 ! 4. Subroutines used :
3097 !
3098 ! Name Type Module Description
3099 ! ----------------------------------------------------------------
3100 ! STRACE Subr. W3SERVMD Subroutine tracing.
3101 ! ----------------------------------------------------------------
3102 !
3103 ! 5. Called by :
3104 !
3105 ! Name Type Module Description
3106 ! ----------------------------------------------------------------
3107 ! ----------------------------------------------------------------
3108 !
3109 ! 6. Error messages :
3110 ! 7. Remarks
3111 ! 8. Structure :
3112 ! 9. Switches :
3113 !
3114 ! !/S Enable subroutine tracing.
3115 !
3116 ! 10. Source code :
3117 !
3118 !/ ------------------------------------------------------------------- /
3119 #ifdef W3_S
3120  USE w3servmd, ONLY: strace
3121 #endif
3122 !
3123  USE w3gdatmd, ONLY: nsea, nseal
3124  USE w3gdatmd, ONLY: gtype, ungtype
3125  USE w3odatmd, only : iaproc, naproc, ntproc
3126  IMPLICIT NONE
3127 !/
3128 !/ ------------------------------------------------------------------- /
3129 !/ Parameter list
3130 !/
3131 !/ ------------------------------------------------------------------- /
3132 !/ Local PARAMETERs
3133 !/
3134 #ifdef W3_S
3135  INTEGER, SAVE :: IENT = 0
3136 #endif
3137 !/
3138 !/ ------------------------------------------------------------------- /
3139 !/
3140  INTEGER ISEA, JSEA
3141  REAL(rkind), allocatable :: ZETA_WORK(:)
3142 #ifdef W3_S
3143  CALL strace (ient, 'VA_SETUP_IOBPD')
3144 #endif
3145 #ifdef W3_DEBUGSTP
3146  WRITE(740+iaproc,*) 'NAPROC=', naproc
3147  WRITE(740+iaproc,*) 'NTPROC=', ntproc
3148  FLUSH(740+iaproc)
3149 #endif
3150  IF (iaproc .le. naproc) THEN
3151 #ifdef W3_DEBUGSTP
3152  WRITE(740+iaproc,*) 'Begin WAVE_SETUP_COMPUTATION'
3153  FLUSH(740+iaproc)
3154 #endif
3155  IF (do_wave_setup) THEN
3156  IF (gtype .EQ. ungtype) THEN
3157  CALL trig_wave_setup_computation
3158  ELSE
3159  CALL fd_wave_setup_computation
3160  END IF
3161  END IF
3162  END IF
3163 #ifdef W3_DEBUGSTP
3164  WRITE(740+iaproc,*) 'End WAVE_SETUP_COMPUTATION'
3165  FLUSH(740+iaproc)
3166 #endif

References do_wave_setup, fd_wave_setup_computation(), w3gdatmd::gtype, w3odatmd::iaproc, w3odatmd::naproc, w3gdatmd::nsea, w3gdatmd::nseal, w3odatmd::ntproc, w3servmd::strace(), trig_wave_setup_computation(), and w3gdatmd::ungtype.

Referenced by w3wavemd::w3wave().

Variable Documentation

◆ do_wave_setup

logical w3wavset::do_wave_setup = .TRUE.

Definition at line 91 of file w3wavset.F90.

91  LOGICAL :: DO_WAVE_SETUP = .true.

Referenced by wave_setup_computation().

◆ ient

integer, save w3wavset::ient = 0

Definition at line 85 of file w3wavset.F90.

85  INTEGER, SAVE :: IENT = 0
w3gdatmd::nseal
integer, pointer nseal
Definition: w3gdatmd.F90:1097
w3gdatmd::iaa
integer, dimension(:), pointer iaa
Definition: w3gdatmd.F90:1124
yownodepool::pdlib_ia
integer, dimension(:), allocatable, target, public pdlib_ia
Definition: yownodepool.F90:81
yowexchangemodule::pdlib_exchange1dreal
subroutine, public pdlib_exchange1dreal(U)
exchange values in U.
Definition: yowexchangeModule.F90:251
w3gdatmd::ygrd
double precision, dimension(:,:), pointer ygrd
Definition: w3gdatmd.F90:1205
yowelementpool
Definition: yowelementpool.F90:38
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3gdatmd::ungtype
integer, parameter ungtype
Definition: w3gdatmd.F90:626
w3wdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: w3wdatmd.F90:18
w3adatmd::cg
real, dimension(:,:), pointer cg
Definition: w3adatmd.F90:575
w3gdatmd::neigh
integer, dimension(:,:), pointer neigh
Definition: w3gdatmd.F90:1105
yownodepool::pdlib_i_diag
integer, dimension(:), allocatable, target, public pdlib_i_diag
Definition: yownodepool.F90:85
w3gdatmd::edges
integer, dimension(:,:), pointer edges
Definition: w3gdatmd.F90:1105
w3adatmd::dw
real, dimension(:), pointer dw
Definition: w3adatmd.F90:584
w3odatmd::ntproc
integer, pointer ntproc
Definition: w3odatmd.F90:457
w3gdatmd::xgrd
double precision, dimension(:,:), pointer xgrd
Definition: w3gdatmd.F90:1205
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3gdatmd::grids
type(grid), dimension(:), allocatable, target grids
Definition: w3gdatmd.F90:1088
w3gdatmd::ny
integer, pointer ny
Definition: w3gdatmd.F90:1097
yownodepool::iplg
integer, dimension(:), allocatable, public iplg
Node local to global mapping.
Definition: yownodepool.F90:116
w3gdatmd::nbedge
integer, pointer nbedge
Definition: w3gdatmd.F90:1104
yownodepool::npa
integer, public npa
number of ghost + resident nodes this partition holds
Definition: yownodepool.F90:99
w3gdatmd::crit_dep_stp
real(8), pointer crit_dep_stp
Definition: w3gdatmd.F90:1409
w3gdatmd::solverthr_stp
real(8), pointer solverthr_stp
Definition: w3gdatmd.F90:1408
yownodepool::pdlib_si
real(rkind), dimension(:), allocatable, target, public pdlib_si
Definition: yownodepool.F90:80
yowelementpool::ne
integer, public ne
number of local elements
Definition: yowelementpool.F90:48
w3gdatmd::mapfs
integer, dimension(:,:), pointer mapfs
Definition: w3gdatmd.F90:1163
yownodepool::y
real(rkind), dimension(:), allocatable, target, public y
Definition: yownodepool.F90:79
yownodepool::pdlib_ja_ie
integer, dimension(:,:,:), allocatable, target, public pdlib_ja_ie
Definition: yownodepool.F90:82
yownodepool
Has data that belong to nodes.
Definition: yownodepool.F90:39
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
yownodepool::x
real(rkind), dimension(:), allocatable, target, public x
coordinates of the local + ghost nodes.
Definition: yownodepool.F90:79
yownodepool::pdlib_ja
integer, dimension(:), allocatable, target, public pdlib_ja
Definition: yownodepool.F90:81
w3parall::synchronize_global_array
subroutine synchronize_global_array(TheVar)
Sync global array in context of PDLIB.
Definition: w3parall.F90:1517
yowdatapool::rtype
integer, save rtype
Definition: yowdatapool.F90:76
w3adatmd::sxy
real, dimension(:), pointer sxy
Definition: w3adatmd.F90:607
w3odatmd
Definition: w3odatmd.F90:3
constants::dwat
real, parameter dwat
DWAT Density of water (kg/m3).
Definition: constants.F90:62
w3gdatmd::mapsf
integer, dimension(:,:), pointer mapsf
Definition: w3gdatmd.F90:1163
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
yownodepool::np
integer, public np
number of nodes, local
Definition: yownodepool.F90:93
yowdatapool::istatus
integer, dimension(mpi_status_size) istatus
MPI Real Type Shpuld be MPI_REAL8.
Definition: yowdatapool.F90:74
yowdatapool
Has fancy data.
Definition: yowdatapool.F90:39
yowexchangemodule
Has only the ghost nodes assign to a neighbor domain.
Definition: yowexchangeModule.F90:39
w3gdatmd::nnz
integer, pointer nnz
Definition: w3gdatmd.F90:1109
w3adatmd::wn
real, dimension(:,:), pointer wn
Definition: w3adatmd.F90:575
w3gdatmd::jaa
integer, dimension(:), pointer jaa
Definition: w3gdatmd.F90:1124
yowelementpool::ine
integer, dimension(:,:), allocatable, target, public ine
number of elements of the augmented domain
Definition: yowelementpool.F90:56
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3gdatmd::gtype
integer, pointer gtype
Definition: w3gdatmd.F90:1094
w3wdatmd::zeta_setup
real, dimension(:), pointer zeta_setup
Definition: w3wdatmd.F90:187
yownodepool::pdlib_tria
real(rkind), dimension(:), allocatable, target, public pdlib_tria
Definition: yownodepool.F90:80
w3gdatmd::si
real(8), dimension(:), pointer si
Definition: w3gdatmd.F90:1122
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
yownodepool::pdlib_ien
real(rkind), dimension(:,:), allocatable, target, public pdlib_ien
Definition: yownodepool.F90:80
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3parall
Parallel routines for implicit solver.
Definition: w3parall.F90:22
w3gdatmd::mapsta
integer, dimension(:,:), pointer mapsta
Definition: w3gdatmd.F90:1163
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61
yownodepool::pdlib_nnz
integer, public pdlib_nnz
Definition: yownodepool.F90:90
w3parall::init_get_isea
subroutine init_get_isea(ISEA, JSEA)
Set ISEA for all schemes.
Definition: w3parall.F90:1398
w3adatmd::syy
real, dimension(:), pointer syy
Definition: w3adatmd.F90:607
w3adatmd::sxx
real, dimension(:), pointer sxx
Definition: w3adatmd.F90:607
w3adatmd::mpi_comm_wcmp
integer, pointer mpi_comm_wcmp
Definition: w3adatmd.F90:676