WAVEWATCH III  beta 0.0.1
w3parall Module Reference

Parallel routines for implicit solver. More...

Functions/Subroutines

subroutine wav_my_wtime (eTime)
 NA. More...
 
subroutine print_my_time (string)
 Print timings. More...
 
subroutine prop_refraction_pr1 (ISEA, DTG, CAD)
 Compute refraction part in matrix. More...
 
subroutine prop_refraction_pr3 (IP, ISEA, DTG, CAD, DoLimiter)
 Compute refraction part in matrix alternative approach. More...
 
subroutine prop_freq_shift (IP, ISEA, CAS, DMM, DTG)
 Compute frequency shift in matrix. More...
 
subroutine prop_freq_shift_m2 (IP, ISEA, CWNB_M2, DWNI_M2, DTG)
 Compute frequency shift alternative approach. More...
 
subroutine synchronize_ipgl_etc_array (IMOD, IsMulti)
 Sync global local arrays. More...
 
subroutine set_up_nseal_nsealm (NSEALout, NSEALMout)
 Setup NSEAL, NSEALM in context of PDLIB. More...
 
subroutine init_get_jsea_isproc (ISEA, JSEA, ISPROC)
 Set JSEA for all schemes. More...
 
subroutine get_jsea_ibelong (ISEA, JSEA, IBELONG)
 Set belongings of JSEA in context of PDLIB. More...
 
subroutine init_get_isea (ISEA, JSEA)
 Set ISEA for all schemes. More...
 
subroutine synchronize_global_array (TheVar)
 Sync global array in context of PDLIB. More...
 

Variables

integer, save ient = 0
 
integer pdlib_nseal
 
integer pdlib_nsealm
 
integer, dimension(:), allocatable jx_to_jsea
 
integer, dimension(:), allocatable isea_to_jsea
 
integer, dimension(:), allocatable listispnextdir
 
integer, dimension(:), allocatable listispprevdir
 
integer, dimension(:), allocatable listispnextfreq
 
integer, dimension(:), allocatable listispprevfreq
 
logical, parameter lsloc = .true.
 
integer, parameter imem = 1
 
real, parameter onesixth = 1.0d0/6.0d0
 
real, parameter onethird = 1.0d0/3.0d0
 
real, parameter zero = 0.0d0
 
real *8, parameter thr8 = TINY(1.d0)
 
real, parameter thr = TINY(1.0)
 

Detailed Description

Parallel routines for implicit solver.

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

Function/Subroutine Documentation

◆ get_jsea_ibelong()

subroutine w3parall::get_jsea_ibelong ( integer, intent(in)  ISEA,
integer, intent(out)  JSEA,
integer, intent(out)  IBELONG 
)

Set belongings of JSEA in context of PDLIB.

Parameters
[in]ISEA
[out]JSEA
[out]IBELONG
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 1271 of file w3parall.F90.

1271  !/ ------------------------------------------------------------------- /
1272  !/
1273  !/ +-----------------------------------+
1274  !/ | WAVEWATCH III NOAA/NCEP |
1275  !/ | |
1276  !/ | Aron Roland (BGS IT&E GmbH) |
1277  !/ | Mathieu Dutour-Sikiric (IRB) |
1278  !/ | |
1279  !/ | FORTRAN 90 |
1280  !/ | Last update : 01-June-2018 |
1281  !/ +-----------------------------------+
1282  !/
1283  !/ 01-June-2018 : Origination. ( version 6.04 )
1284  !/
1285  ! 1. Purpose : Set belongings of jsea in context of pdlib
1286  ! 2. Method :
1287  ! 3. Parameters :
1288  !
1289  ! Parameter list
1290  ! ----------------------------------------------------------------
1291  ! ----------------------------------------------------------------
1292  !
1293  ! 4. Subroutines used :
1294  !
1295  ! Name Type Module Description
1296  ! ----------------------------------------------------------------
1297  ! STRACE Subr. W3SERVMD Subroutine tracing.
1298  ! ----------------------------------------------------------------
1299  !
1300  ! 5. Called by :
1301  !
1302  ! Name Type Module Description
1303  ! ----------------------------------------------------------------
1304  ! ----------------------------------------------------------------
1305  !
1306  ! 6. Error messages :
1307  ! 7. Remarks
1308  ! 8. Structure :
1309  ! 9. Switches :
1310  !
1311  ! !/S Enable subroutine tracing.
1312  !
1313  ! 10. Source code :
1314  !
1315  !/ ------------------------------------------------------------------- /
1316 #ifdef W3_S
1317  USE w3servmd, ONLY: strace
1318 #endif
1319  !/
1320  USE w3odatmd, ONLY: outpts, iaproc, naproc
1321  USE w3gdatmd, ONLY: gtype, ungtype, mapsf
1322  USE constants, ONLY : lpdlib
1323 #ifdef W3_PDLIB
1325  use yownodepool, only: ipgl, iplg
1326 #endif
1327  IMPLICIT NONE
1328  !/ ------------------------------------------------------------------- /
1329  !/ Parameter list
1330  !/
1331  !/ ------------------------------------------------------------------- /
1332  !/ Local parameters
1333  !/
1334 #ifdef W3_S
1335  INTEGER, SAVE :: IENT = 0
1336 #endif
1337  !/
1338  !/ ------------------------------------------------------------------- /
1339  !/
1340  INTEGER, intent(in) :: ISEA
1341  INTEGER, intent(out) :: JSEA, IBELONG
1342  INTEGER ISPROC, IX, JX
1343 #ifdef W3_S
1344  CALL strace (ient, 'GET_JSEA_IBELONG')
1345 #endif
1346  IF (.NOT. lpdlib) THEN
1347  jsea = 1 + (isea-1)/naproc
1348  isproc = isea - (jsea-1)*naproc
1349  IF (isproc .eq. iaproc) THEN
1350  ibelong=1
1351  ELSE
1352  ibelong=0
1353  END IF
1354  ELSE
1355 #ifdef W3_PDLIB
1356  IF (gtype .ne. ungtype) THEN
1357  jsea = 1 + (isea-1)/naproc
1358  isproc = isea - (jsea-1)*naproc
1359  IF (isproc .eq. iaproc) THEN
1360  ibelong=1
1361  ELSE
1362  ibelong=0
1363  END IF
1364  ELSE
1365  IF (iaproc .le. naproc) THEN
1366  ix = mapsf(isea,1)
1367  jx = ipgl_npa(ix)
1368  IF (jx .eq. 0) THEN
1369  ibelong=0
1370  jsea=-1
1371  ELSE
1372  ibelong=1
1373  jsea = jx_to_jsea(jx)
1374  END IF
1375  ELSE
1376  ibelong=0
1377  jsea=-1
1378  END IF
1379  ENDIF
1380 #endif
1381  ENDIF
1382  !/
1383  !/ End of INIT_GET_ISEA ---------------------------------------------- /
1384  !/

References w3gdatmd::gtype, w3odatmd::iaproc, yownodepool::ipgl, yowrankmodule::ipgl_npa, yowrankmodule::ipgl_to_proc, yowrankmodule::ipgl_tot, yownodepool::iplg, jx_to_jsea, constants::lpdlib, w3gdatmd::mapsf, w3odatmd::naproc, w3odatmd::outpts, w3servmd::strace(), and w3gdatmd::ungtype.

Referenced by pdlib_field_vec::unst_pdlib_read_from_file(), pdlib_field_vec::unst_pdlib_write_to_file(), w3updtmd::w3uini(), and w3updtmd::w3ulev().

◆ init_get_isea()

subroutine w3parall::init_get_isea ( integer, intent(out)  ISEA,
integer, intent(in)  JSEA 
)

Set ISEA for all schemes.

Parameters
[out]ISEA
[in]JSEA
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 1398 of file w3parall.F90.

1398  !/ ------------------------------------------------------------------- /
1399  !/
1400  !/ +-----------------------------------+
1401  !/ | WAVEWATCH III NOAA/NCEP |
1402  !/ | |
1403  !/ | Aron Roland (BGS IT&E GmbH) |
1404  !/ | Mathieu Dutour-Sikiric (IRB) |
1405  !/ | |
1406  !/ | FORTRAN 90 |
1407  !/ | Last update : 01-June-2018 |
1408  !/ +-----------------------------------+
1409  !/
1410  !/ 01-June-2018 : Origination. ( version 6.04 )
1411  !/
1412  ! 1. Purpose : Set Isea for all schemes
1413  ! 2. Method :
1414  ! 3. Parameters :
1415  !
1416  ! Parameter list
1417  ! ----------------------------------------------------------------
1418  ! ----------------------------------------------------------------
1419  !
1420  ! 4. Subroutines used :
1421  !
1422  ! Name Type Module Description
1423  ! ----------------------------------------------------------------
1424  ! STRACE Subr. W3SERVMD Subroutine tracing.
1425  ! ----------------------------------------------------------------
1426  !
1427  ! 5. Called by :
1428  !
1429  ! Name Type Module Description
1430  ! ----------------------------------------------------------------
1431  ! ----------------------------------------------------------------
1432  !
1433  ! 6. Error messages :
1434  ! 7. Remarks
1435  ! 8. Structure :
1436  ! 9. Switches :
1437  !
1438  ! !/S Enable subroutine tracing.
1439  !
1440  ! 10. Source code :
1441  !
1442  !/ ------------------------------------------------------------------- /
1443 #ifdef W3_S
1444  USE w3servmd, ONLY: strace
1445 #endif
1446  !/
1447  USE w3odatmd, ONLY: outpts, iaproc, naproc
1448  USE w3gdatmd, ONLY: gtype, ungtype
1449  USE constants, ONLY : lpdlib
1450 #ifdef W3_PDLIB
1451  USE yownodepool, ONLY: iplg
1452 #endif
1453  !/ ------------------------------------------------------------------- /
1454  !/ Parameter list
1455  !/
1456  !/ ------------------------------------------------------------------- /
1457  !/ Local parameters
1458  !/
1459  !/ ------------------------------------------------------------------- /
1460  !/
1461  !/
1462  !/ ------------------------------------------------------------------- /
1463  !/
1464  !/ ------------------------------------------------------------------- /
1465  !
1466  USE w3odatmd, ONLY: outpts, iaproc, naproc
1467  USE w3gdatmd, ONLY: gtype, ungtype
1468  USE constants, ONLY : lpdlib
1469 #ifdef W3_PDLIB
1470  USE yownodepool, ONLY: iplg
1471 #endif
1472  IMPLICIT NONE
1473 #ifdef W3_S
1474  INTEGER, SAVE :: IENT = 0
1475 #endif
1476  INTEGER, intent(in) :: JSEA
1477  INTEGER, intent(out) :: ISEA
1478 #ifdef W3_S
1479  CALL strace (ient, 'INIT_GET_ISEA')
1480 #endif
1481 #ifdef W3_SHRD
1482  isea = jsea
1483 #endif
1484 #ifdef W3_DIST
1485  IF (.NOT. lpdlib) THEN
1486  isea = iaproc + (jsea-1)*naproc
1487  ELSE
1488 #endif
1489 #ifdef W3_PDLIB
1490  IF (gtype .eq. ungtype) THEN
1491  isea = iplg(jsea)
1492  ELSE
1493  isea = iaproc + (jsea-1)*naproc
1494  ENDIF
1495 #endif
1496 #ifdef W3_DIST
1497  ENDIF
1498 #endif
1499  !/
1500  !/ End of INIT_GET_ISEA ------------------------------------------------ /
1501  !/

References w3gdatmd::gtype, w3odatmd::iaproc, yownodepool::iplg, constants::lpdlib, w3odatmd::naproc, w3odatmd::outpts, w3servmd::strace(), and w3gdatmd::ungtype.

Referenced by pdlib_w3profsmd::block_solver_init(), w3iogomd::calc_u3stokes(), w3iogomd::calc_wbt(), pdlib_w3profsmd::check_array_integral_nx_r8_maxfunct(), w3oacpmd::cpl_oasis_define(), w3wavset::differentiate_xydir_mapsta(), pdlib_field_vec::do_output_exchanges(), pdlib_w3profsmd::scal_integral_print_general(), pdlib_w3profsmd::set_iobdp_pdlib(), pdlib_w3profsmd::set_iobpa_pdlib(), w3iogomd::skewness(), w3wavset::trig_compute_lh_stress(), w3wavset::trig_wave_setup_computation(), pdlib_field_vec::unst_pdlib_write_to_file(), w3iosfmd::w3cprt(), w3wavemd::w3gath(), w3initmd::w3init(), w3iorsmd::w3iors(), w3iogomd::w3outg(), w3wavemd::w3scat(), w3snl5md::w3snl5(), w3updtmd::w3uice(), w3updtmd::w3uini(), w3updtmd::w3ulev(), w3wavemd::w3wave(), wminiomd::wmiohg(), wminiomd::wmiohs(), and pdlib_w3profsmd::write_var_to_text_file().

◆ init_get_jsea_isproc()

subroutine w3parall::init_get_jsea_isproc ( integer, intent(in)  ISEA,
integer, intent(out)  JSEA,
integer, intent(out)  ISPROC 
)

Set JSEA for all schemes.

Parameters
[in]ISEA
[out]JSEA
[out]ISPROC
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 1163 of file w3parall.F90.

1163  !/ ------------------------------------------------------------------- /
1164  !/
1165  !/ +-----------------------------------+
1166  !/ | WAVEWATCH III NOAA/NCEP |
1167  !/ | |
1168  !/ | Aron Roland (BGS IT&E GmbH) |
1169  !/ | Mathieu Dutour-Sikiric (IRB) |
1170  !/ | |
1171  !/ | FORTRAN 90 |
1172  !/ | Last update : 01-June-2018 |
1173  !/ +-----------------------------------+
1174  !/
1175  !/ 01-June-2018 : Origination. ( version 6.04 )
1176  !/
1177  ! 1. Purpose : Set Jsea for all schemes
1178  ! 2. Method :
1179  ! 3. Parameters :
1180  !
1181  ! Parameter list
1182  ! ----------------------------------------------------------------
1183  ! ----------------------------------------------------------------
1184  !
1185  ! 4. Subroutines used :
1186  !
1187  ! Name Type Module Description
1188  ! ----------------------------------------------------------------
1189  ! STRACE Subr. W3SERVMD Subroutine tracing.
1190  ! ----------------------------------------------------------------
1191  !
1192  ! 5. Called by :
1193  !
1194  ! Name Type Module Description
1195  ! ----------------------------------------------------------------
1196  ! ----------------------------------------------------------------
1197  !
1198  ! 6. Error messages :
1199  ! 7. Remarks
1200  ! 8. Structure :
1201  ! 9. Switches :
1202  !
1203  ! !/S Enable subroutine tracing.
1204  !
1205  ! 10. Source code :
1206  !
1207  !/ ------------------------------------------------------------------- /
1208 #ifdef W3_S
1209  USE w3servmd, ONLY: strace
1210 #endif
1211  !/
1212  USE w3odatmd, ONLY: outpts, iaproc, naproc
1213  USE w3gdatmd, ONLY: gtype, ungtype, mapsf
1214  USE constants, ONLY : lpdlib
1215 #ifdef W3_PDLIB
1216  USE yowrankmodule, only : ipgl_to_proc, ipgl_tot
1217  use yownodepool, only: ipgl, iplg
1218 #endif
1219  IMPLICIT NONE
1220  !/ ------------------------------------------------------------------- /
1221  !/ Parameter list
1222  !/
1223  !/ ------------------------------------------------------------------- /
1224  !/ Local parameters
1225  !/
1226 #ifdef W3_S
1227  INTEGER, SAVE :: IENT = 0
1228 #endif
1229  !/
1230  !/ ------------------------------------------------------------------- /
1231  INTEGER, intent(in) :: ISEA
1232  INTEGER, intent(out) :: JSEA, ISPROC
1233  INTEGER IP_glob
1234 #ifdef W3_S
1235  CALL strace (ient, 'INIT_GET_JSEA_ISPROC')
1236 #endif
1237 
1238 #ifdef W3_PDLIB
1239  IF (.NOT. lpdlib) THEN
1240 #endif
1241  jsea = 1 + (isea-1)/naproc
1242  isproc = isea - (jsea-1)*naproc
1243 #ifdef W3_PDLIB
1244  ELSE
1245  ip_glob = mapsf(isea,1)
1246  IF (iaproc .le. naproc) THEN
1247  jsea = isea_to_jsea(isea)
1248  ELSE
1249  jsea = -1
1250  END IF
1251  isproc = ipgl_to_proc(ip_glob)
1252  ENDIF
1253 #endif
1254  !/
1255  !/ End of JACOBI_INIT ------------------------------------------------ /
1256  !/

References w3gdatmd::gtype, w3odatmd::iaproc, yownodepool::ipgl, yowrankmodule::ipgl_to_proc, yowrankmodule::ipgl_tot, yownodepool::iplg, isea_to_jsea, constants::lpdlib, w3gdatmd::mapsf, w3odatmd::naproc, w3odatmd::outpts, w3servmd::strace(), and w3gdatmd::ungtype.

Referenced by pdlib_field_vec::do_output_exchanges(), pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block(), pdlib_w3profsmd::pdlib_w3xypfsfct2(), pdlib_w3profsmd::pdlib_w3xypfsn2(), pdlib_w3profsmd::pdlib_w3xypfspsi2(), w3iosfmd::w3cprt(), w3initmd::w3init(), w3iorsmd::w3iors(), w3iosfmd::w3iosf(), w3iotrmd::w3iotr(), w3initmd::w3mpio(), w3initmd::w3mpip(), w3wavemd::w3nmin(), w3updtmd::w3uice(), w3updtmd::w3uini(), w3updtmd::w3ulev(), wmgridmd::wmghgh(), wmgridmd::wmglow(), and wminiomd::wmiobs().

◆ print_my_time()

subroutine w3parall::print_my_time ( character(*), intent(in)  string)

Print timings.

Parameters
[in]string
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 200 of file w3parall.F90.

200  !/
201  !/ +-----------------------------------+
202  !/ | WAVEWATCH III NOAA/NCEP |
203  !/ | |
204  !/ | Aron Roland (BGS IT&E GmbH) |
205  !/ | Mathieu Dutour-Sikiric (IRB) |
206  !/ | |
207  !/ | FORTRAN 90 |
208  !/ | Last update : 01-June-2018 |
209  !/ +-----------------------------------+
210  !/
211  !/ 01-June-2018 : Origination. ( version 6.04 )
212  !/
213  ! 1. Purpose : Print timings
214  ! 2. Method :
215  ! 3. Parameters :
216  !
217  ! Parameter list
218  ! ----------------------------------------------------------------
219  ! ----------------------------------------------------------------
220  !
221  ! 4. Subroutines used :
222  !
223  ! Name Type Module Description
224  ! ----------------------------------------------------------------
225  ! STRACE Subr. W3SERVMD Subroutine tracing.
226  ! ----------------------------------------------------------------
227  !
228  ! 5. Called by :
229  !
230  ! Name Type Module Description
231  ! ----------------------------------------------------------------
232  ! ----------------------------------------------------------------
233  !
234  ! 6. Error messages :
235  ! 7. Remarks
236  ! 8. Structure :
237  ! 9. Switches :
238  !
239  ! !/S Enable subroutine tracing.
240  !
241  ! 10. Source code :
242  !
243  !/ ------------------------------------------------------------------- /
244 #ifdef W3_S
245  USE w3servmd, ONLY: strace
246 #endif
247  USE w3odatmd, ONLY : iaproc
248  IMPLICIT NONE
249  !/
250  !/ ------------------------------------------------------------------- /
251  !/ Parameter list
252  !/
253  !/ ------------------------------------------------------------------- /
254  !/ Local parameters
255  !/
256 #ifdef W3_S
257  INTEGER, SAVE :: IENT = 0
258 #endif
259  !/
260  !/ ------------------------------------------------------------------- /
261  !
262  character(*), intent(in) :: string
263  REAL(8) :: eTime
264 #ifdef W3_S
265  CALL strace (ient, 'PRINT_MY_TIME')
266 #endif
267  CALL wav_my_wtime(etime)
268  WRITE(740+iaproc,*) 'TIMING time=', etime, ' at step ', string
269  !/
270  !/ End of JACOBI_INIT ------------------------------------------------ /
271  !/

References w3odatmd::iaproc, w3servmd::strace(), and wav_my_wtime().

Referenced by pdlib_field_vec::unst_pdlib_read_from_file(), w3initmd::w3init(), w3iorsmd::w3iors(), and w3wavemd::w3wave().

◆ prop_freq_shift()

subroutine w3parall::prop_freq_shift ( integer, intent(in)  IP,
integer, intent(in)  ISEA,
real, dimension(nspec), intent(out)  CAS,
real, dimension(0:nk2), intent(out)  DMM,
real, intent(in)  DTG 
)

Compute frequency shift in matrix.

Parameters
[in]IP
[in]ISEA
[out]CAS
[out]DMM
[in]DTG
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 609 of file w3parall.F90.

609  !/
610  !/ +-----------------------------------+
611  !/ | WAVEWATCH III NOAA/NCEP |
612  !/ | |
613  !/ | Aron Roland (BGS IT&E GmbH) |
614  !/ | Mathieu Dutour-Sikiric (IRB) |
615  !/ | |
616  !/ | FORTRAN 90 |
617  !/ | Last update : 01-June-2018 |
618  !/ +-----------------------------------+
619  !/
620  !/ 01-June-2018 : Origination. ( version 6.04 )
621  !/
622  ! 1. Purpose : Compute freq. shift in matrix
623  ! 2. Method :
624  ! 3. Parameters :
625  !
626  ! Parameter list
627  ! ----------------------------------------------------------------
628  ! ----------------------------------------------------------------
629  !
630  ! 4. Subroutines used :
631  !
632  ! Name Type Module Description
633  ! ----------------------------------------------------------------
634  ! STRACE Subr. W3SERVMD Subroutine tracing.
635  ! ----------------------------------------------------------------
636  !
637  ! 5. Called by :
638  !
639  ! Name Type Module Description
640  ! ----------------------------------------------------------------
641  ! ----------------------------------------------------------------
642  !
643  ! 6. Error messages :
644  ! 7. Remarks
645  ! 8. Structure :
646  ! 9. Switches :
647  !
648  ! !/S Enable subroutine tracing.
649  !
650  ! 10. Source code :
651  !
652  !/ ------------------------------------------------------------------- /
653 #ifdef W3_S
654  USE w3servmd, ONLY: strace
655 #endif
656  USE constants, ONLY : lpdlib
657  USE w3gdatmd, ONLY: nk, nk2, nth, nspec, sig, dsip, ecos, esin, &
658  ec2, esc, es2, fachfa, mapwn, flcth, flck, &
659  ctmax, dmin, dth, mapsf
660  USE w3adatmd, ONLY: cg, wn, dcxdx, dcxdy, dcydx, dcydy, cx, cy, dddx, dddy, dw
661  USE w3odatmd, only : iaproc
662  IMPLICIT NONE
663  !/ Parameter list
664  !/
665  !/ ------------------------------------------------------------------- /
666  !/ Local parameters
667  !/
668 #ifdef W3_S
669  INTEGER, SAVE :: IENT = 0
670 #endif
671  INTEGER, intent(in) :: ISEA, IP
672  REAL, intent(out) :: DMM(0:NK2)
673  REAL, intent(in) :: DTG
674  REAL, intent(out) :: CAS(NSPEC)
675  REAL :: DB(NK2), DSDD(0:NK+1)
676  REAL :: eDCXDX, eDCXDY, eDCYDX, eDCYDY, eCX, eCY, eDDDX, EDDDY
677  REAL :: DCXX, DCXYYX, DCYY, FKD, FACK
678  REAL :: VELNOFILT, VELFAC, DEPTH
679  REAL :: CFLK(NK2,NTH), FKC(NTH), FKD0
680  INTEGER :: IK, ITH, ISP, IY, IX
681 #ifdef W3_S
682  CALL strace (ient, 'PROP_FREQ_SHIFT')
683 #endif
684  !
685  IF (lpdlib) THEN
686  edcxdx = dcxdx(1,ip)
687  edcxdy = dcxdy(1,ip)
688  edcydx = dcydx(1,ip)
689  edcydy = dcydy(1,ip)
690  edddx = dddx(1,ip)
691  edddy = dddy(1,ip)
692  ELSE
693  ix=mapsf(isea,1)
694  iy=mapsf(isea,2)
695  edcxdx=dcxdx(iy,ix)
696  edcxdy=dcxdy(iy,ix)
697  edcydx=dcydx(iy,ix)
698  edcydy=dcydy(iy,ix)
699  edddx=dddx(iy,ix)
700  edddy=dddy(iy,ix)
701  ENDIF
702  ecx=cx(isea)
703  ecy=cy(isea)
704  dcxx = - edcxdx
705  dcxyyx = - ( edcxdy + edcydx )
706  dcyy = - edcydy
707  fkd = ( ecx*edddx + ecy*edddy )
708  fack = dtg
709  DO ith=1, nth
710  fkc(ith) = ec2(ith)*dcxx + esc(ith)*dcxyyx + es2(ith)*dcyy
711  END DO
712  DO ik=0, nk
713  db(ik+1) = dsip(ik) / cg(ik,isea)
714  dmm(ik+1) = dble(wn(ik+1,isea) - wn(ik,isea))
715  END DO
716  db(nk+2) = dsip(nk+1) / cg(nk+1,isea)
717  dmm(nk+2) = zero
718  dmm(0)=dmm(1)
719  !
720  depth = max( dmin , dw(isea) )
721  DO ik=0, nk+1
722  IF ( depth*wn(ik,isea) .LT. 5. ) THEN
723  dsdd(ik) = max( 0. , cg(ik,isea)*wn(ik,isea)-0.5*sig(ik) ) / depth
724  ELSE
725  dsdd(ik) = 0.
726  END IF
727  END DO
728  DO ik=0, nk+1
729  fkd0 = fkd / cg(ik,isea) * dsdd(ik)
730  velfac = fack/db(ik+1)
731  DO ith=1, nth
732  velnofilt = ( fkd0 + wn(ik,isea)*fkc(ith) ) * velfac
733  cflk(ik+1,ith) = velnofilt/velfac
734  END DO
735  END DO
736  DO ik=1,nk
737  DO ith=1,nth
738  isp=ith + (ik-1)*nth
739  cas(isp)=dble(cflk(ik,ith))
740  END DO
741  END DO
742  !/
743  !/ End of JACOBI_INIT ------------------------------------------------ /
744  !/

References w3adatmd::cg, w3gdatmd::ctmax, w3adatmd::cx, w3adatmd::cy, w3adatmd::dcxdx, w3adatmd::dcxdy, w3adatmd::dcydx, w3adatmd::dcydy, w3adatmd::dddx, w3adatmd::dddy, w3gdatmd::dmin, w3gdatmd::dsip, w3gdatmd::dth, w3adatmd::dw, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::es2, w3gdatmd::esc, w3gdatmd::esin, w3gdatmd::fachfa, w3gdatmd::flck, w3gdatmd::flcth, w3odatmd::iaproc, constants::lpdlib, w3gdatmd::mapsf, w3gdatmd::mapwn, w3gdatmd::nk, w3gdatmd::nk2, w3gdatmd::nspec, w3gdatmd::nth, w3gdatmd::sig, w3servmd::strace(), w3adatmd::wn, and zero.

Referenced by pdlib_w3profsmd::calcarray_jacobi_spectral_1(), and pdlib_w3profsmd::calcarray_jacobi_spectral_2().

◆ prop_freq_shift_m2()

subroutine w3parall::prop_freq_shift_m2 ( integer, intent(in)  IP,
integer, intent(in)  ISEA,
real, dimension(1-nth:nspec), intent(out)  CWNB_M2,
real, dimension(nk), intent(out)  DWNI_M2,
real, intent(in)  DTG 
)

Compute frequency shift alternative approach.

Parameters
[in]IP
[in]ISEA
[out]CWNB_M2
[out]DWNI_M2
[in]DTG
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 761 of file w3parall.F90.

761  !/
762  !/ +-----------------------------------+
763  !/ | WAVEWATCH III NOAA/NCEP |
764  !/ | |
765  !/ | Aron Roland (BGS IT&E GmbH) |
766  !/ | Mathieu Dutour-Sikiric (IRB) |
767  !/ | |
768  !/ | FORTRAN 90 |
769  !/ | Last update : 01-June-2018 |
770  !/ +-----------------------------------+
771  !/
772  !/ 01-June-2018 : Origination. ( version 6.04 )
773  !/
774  ! 1. Purpose : Compute freq. shift alternative approach
775  ! 2. Method :
776  ! 3. Parameters :
777  !
778  ! Parameter list
779  ! ----------------------------------------------------------------
780  ! ----------------------------------------------------------------
781  !
782  ! 4. Subroutines used :
783  !
784  ! Name Type Module Description
785  ! ----------------------------------------------------------------
786  ! STRACE Subr. W3SERVMD Subroutine tracing.
787  ! ----------------------------------------------------------------
788  !
789  ! 5. Called by :
790  !
791  ! Name Type Module Description
792  ! ----------------------------------------------------------------
793  ! ----------------------------------------------------------------
794  !
795  ! 6. Error messages :
796  ! 7. Remarks
797  ! 8. Structure :
798  ! 9. Switches :
799  !
800  ! !/S Enable subroutine tracing.
801  !
802  ! 10. Source code :
803  !
804  !/ ------------------------------------------------------------------- /
805  !
806 #ifdef W3_S
807  USE w3servmd, ONLY: strace
808 #endif
809  USE constants, ONLY : lpdlib
810  USE w3gdatmd, ONLY: nk, nk2, nth, nspec, sig, dsip, ecos, esin, &
811  ec2, esc, es2, fachfa, mapwn, flcth, flck, &
812  ctmax, dmin, dth, mapsf
813  USE w3adatmd, ONLY: cg, wn, dcxdx, dcxdy, dcydx, dcydy, cx, cy, dddx, dddy, dw
814  USE w3odatmd, only : iaproc
815 
816  IMPLICIT NONE
817 
818  !/ ------------------------------------------------------------------- /
819  !/ Local parameters
820  !/
821 #ifdef W3_S
822  INTEGER, SAVE :: IENT = 0
823 #endif
824 
825  INTEGER, intent(in) :: ISEA, IP
826  REAL, intent(out) :: CWNB_M2(1-NTH:NSPEC)
827  REAL, intent(out) :: DWNI_M2(NK)
828  REAL, intent(in) :: DTG
829  !
830  REAL :: eDCXDX, eDCXDY, eDCYDX, eDCYDY, eCX, eCY, eDDDX, EDDDY
831  REAL :: DCXX, DCXYYX, DCYY, FKD, FACK
832  REAL :: DEPTH
833  REAL :: FKC(NTH), FKD0
834  REAL :: VCWN(1-NTH:NSPEC+NTH)
835  REAL :: DSDD(0:NK+1)
836  REAL :: sumDiff, sumDiff1, sumDiff2, sumDiff3
837  REAL :: sumDiff0, sumDiff4, sumDiff5
838  INTEGER :: IK, ITH, ISP, IY, IX
839 
840  !/ ------------------------------------------------------------------- /
841 #ifdef W3_S
842  CALL strace (ient, 'PROP_FREQ_SHIFT_M2')
843 #endif
844 
845  IF (lpdlib) THEN
846  edcxdx = dcxdx(1,ip)
847  edcxdy = dcxdy(1,ip)
848  edcydx = dcydx(1,ip)
849  edcydy = dcydy(1,ip)
850  edddx = dddx(1,ip)
851  edddy = dddy(1,ip)
852  ELSE
853  ix=mapsf(isea,1)
854  iy=mapsf(isea,2)
855  edcxdx=dcxdx(iy,ix)
856  edcxdy=dcxdy(iy,ix)
857  edcydx=dcydx(iy,ix)
858  edcydy=dcydy(iy,ix)
859  edddx=dddx(iy,ix)
860  edddy=dddy(iy,ix)
861  ENDIF
862 
863  ecx = cx(isea)
864  ecy = cy(isea)
865  fack = dtg
866  dcxx = - fack * edcxdx
867  dcxyyx = - fack * ( edcxdy + edcydx )
868  dcyy = - fack * edcydy
869  fkd = fack * ( ecx*edddx + ecy*edddy )
870 
871  DO ith=1, nth
872  fkc(ith) = ec2(ith)*dcxx + esc(ith)*dcxyyx + es2(ith)*dcyy
873  END DO
874  !
875  depth = max( dmin , dw(isea) )
876  DO ik=0, nk+1
877  IF ( depth*wn(ik,isea) .LT. 5. ) THEN
878  dsdd(ik) = max( 0. , cg(ik,isea)*wn(ik,isea)-0.5*sig(ik) ) / depth
879  ELSE
880  dsdd(ik) = 0.
881  END IF
882  END DO
883  isp = -nth
884  DO ik=0, nk+1
885  fkd0 = fkd / cg(ik,isea) * dsdd(ik)
886  DO ith=1, nth
887  isp = isp + 1
888  vcwn(isp) = fkd0 + wn(ik,isea)*fkc(ith)
889  END DO
890  END DO
891 
892  sumdiff=0
893  DO isp=1-nth,nspec
894  cwnb_m2(isp) = dble(0.5 * ( vcwn(isp) + vcwn(isp+nth) ))
895  sumdiff = sumdiff + max(cwnb_m2(isp), zero)
896  END DO
897  DO ik=1,nk
898  dwni_m2(ik) = dble( cg(ik,isea) / dsip(ik) )
899  END DO
900  !/
901  !/ End of JACOBI_INIT ------------------------------------------------ /
902  !/

References w3adatmd::cg, w3gdatmd::ctmax, w3adatmd::cx, w3adatmd::cy, w3adatmd::dcxdx, w3adatmd::dcxdy, w3adatmd::dcydx, w3adatmd::dcydy, w3adatmd::dddx, w3adatmd::dddy, w3gdatmd::dmin, w3gdatmd::dsip, w3gdatmd::dth, w3adatmd::dw, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::es2, w3gdatmd::esc, w3gdatmd::esin, w3gdatmd::fachfa, w3gdatmd::flck, w3gdatmd::flcth, w3odatmd::iaproc, constants::lpdlib, w3gdatmd::mapsf, w3gdatmd::mapwn, w3gdatmd::nk, w3gdatmd::nk2, w3gdatmd::nspec, w3gdatmd::nth, w3gdatmd::sig, w3servmd::strace(), w3adatmd::wn, and zero.

Referenced by pdlib_w3profsmd::calcarray_jacobi_spectral_1(), and pdlib_w3profsmd::calcarray_jacobi_spectral_2().

◆ prop_refraction_pr1()

subroutine w3parall::prop_refraction_pr1 ( integer, intent(in)  ISEA,
real, intent(in)  DTG,
real, dimension(nspec), intent(out)  CAD 
)

Compute refraction part in matrix.

Parameters
[in]ISEA
[in]DTG
[out]CAD
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 286 of file w3parall.F90.

286  !/
287  !/ +-----------------------------------+
288  !/ | WAVEWATCH III NOAA/NCEP |
289  !/ | |
290  !/ | Aron Roland (BGS IT&E GmbH) |
291  !/ | Mathieu Dutour-Sikiric (IRB) |
292  !/ | |
293  !/ | FORTRAN 90 |
294  !/ | Last update : 01-June-2018 |
295  !/ +-----------------------------------+
296  !/
297  !/ 01-June-2018 : Origination. ( version 6.04 )
298  !/
299  ! 1. Purpose : Compute refraction part in matrix
300  ! 2. Method :
301  ! 3. Parameters :
302  !
303  ! Parameter list
304  ! ----------------------------------------------------------------
305  ! ----------------------------------------------------------------
306  !
307  ! 4. Subroutines used :
308  !
309  ! Name Type Module Description
310  ! ----------------------------------------------------------------
311  ! STRACE Subr. W3SERVMD Subroutine tracing.
312  ! ----------------------------------------------------------------
313  !
314  ! 5. Called by :
315  !
316  ! Name Type Module Description
317  ! ----------------------------------------------------------------
318  ! ----------------------------------------------------------------
319  !
320  ! 6. Error messages :
321  ! 7. Remarks
322  ! 8. Structure :
323  ! 9. Switches :
324  !
325  ! !/S Enable subroutine tracing.
326  !
327  ! 10. Source code :
328  !
329  !/ ------------------------------------------------------------------- /
330 #ifdef W3_S
331  USE w3servmd, ONLY: strace
332 #endif
333  USE w3gdatmd, ONLY: nk, nk2, nth, nspec, sig, dsip, ecos, esin, &
334  ec2, esc, es2, fachfa, mapwn, flcth, flck, &
336  USE w3adatmd, ONLY: cg, wn, dcxdx, dcxdy, dcydx, dcydy, dddx, &
337  dddy, dw
338 #ifdef W3_REFRX
339  USE w3adatmd, ONLY: dcdx, dcdy
340 #endif
341  USE w3idatmd, ONLY: flcur
342  USE w3odatmd, only : iaproc
343  IMPLICIT NONE
344  !/
345  !/ ------------------------------------------------------------------- /
346  !/ Parameter list
347  !/
348  !/ ------------------------------------------------------------------- /
349  !/ Local parameters
350  !/
351 #ifdef W3_S
352  INTEGER, SAVE :: IENT = 0
353 #endif
354  !/
355  !/ ------------------------------------------------------------------- /
356  !/
357  !/
358  REAL, intent(out) :: CAD(NSPEC)
359  INTEGER, intent(in) :: ISEA
360  REAL, intent(in) :: DTG
361  INTEGER :: ISP, IK, ITH, IX, IY
362  REAL :: FRK(NK), FRG(NK), DSDD(0:NK+1)
363  REAL :: FACTH, DCXY, DCYX, DCXXYY, DTTST
364  REAL :: eDCXDX, eDCXDY, eDCYDX, eDCYDY, eDDDX, eDDDY, eCTHG0
365  REAL :: VCFLT(NSPEC), DEPTH, FDG
366  REAL :: FDDMAX
367 #ifdef W3_S
368  CALL strace (ient, 'PROP_REFRACTION_PR1')
369 #endif
370  ix=mapsf(isea,1)
371  iy=mapsf(isea,2)
372  edddx=dddx(iy,ix)
373  edddy=dddy(iy,ix)
374  ecthg0 = cthg0s(isea)
375  facth = dtg / dth
376  !
377  fdg = facth * ecthg0
378  depth = max( dmin , dw(isea) )
379  DO ik=0, nk+1
380  IF ( depth*wn(ik,isea) .LT. 5. ) THEN
381  dsdd(ik) = max( 0. , cg(ik,isea)*wn(ik,isea)-0.5*sig(ik) ) / depth
382  ELSE
383  dsdd(ik) = 0.
384  END IF
385  END DO
386  fddmax=0
387  DO ith=1, nth
388  fddmax = max( fddmax , abs(esin(ith)*edddx - ecos(ith)*edddy ) )
389  END DO
390  DO ik=1, nk
391  frk(ik) = facth * dsdd(ik) / wn(ik,isea)
392  !FRK(IK) = FRK(IK) / MAX ( 1. , FRK(IK)*FDDMAX/CTMAX )
393  frg(ik) = fdg * cg(ik,isea)
394  END DO
395  DO isp=1, nspec
396  vcflt(isp) = frg(mapwn(isp)) * ecos(isp) + &
397  frk(mapwn(isp)) * ( esin(isp)*edddx - ecos(isp)*edddy )
398  END DO
399  !
400 #ifdef W3_REFRX
401  ! 3.c @C/@x refraction and great-circle propagation
402  vcflt = 0.
403  frk = 0.
404  DO ik=1, nk
405  frk(ik) = facth * cg(ik,isea) * wn(ik,isea) / sig(ik)
406  END DO
407  DO isp=1, nspec
408  vcflt(isp) = frg(mapwn(isp)) * ecos(isp) &
409  + frk(mapwn(isp)) * ( esin(isp)*dcdx(isp,1,mapwn(isp)) &
410  - ecos(isp)*dcdy(isp,1,mapwn(isp)) )
411  END DO
412 #endif
413  !
414  IF ( flcur ) THEN
415  edcxdx=dcxdx(iy,ix)
416  edcxdy=dcxdy(iy,ix)
417  edcydx=dcydx(iy,ix)
418  edcydy=dcydy(iy,ix)
419  dcyx = facth * edcydx
420  dcxxyy = facth * ( edcxdx - edcydy )
421  dcxy = facth * edcxdy
422  DO isp=1, nspec
423  vcflt(isp) = vcflt(isp) + es2(isp)*dcyx + esc(isp)*dcxxyy - ec2(isp)*dcxy
424  END DO
425  END IF
426  DO isp=1,nspec
427  cad(isp)=dble(vcflt(isp))
428  END DO
429  !/
430  !/ End of JACOBI_INIT ------------------------------------------------ /
431  !/

References w3adatmd::cg, w3gdatmd::cthg0s, w3gdatmd::ctmax, w3adatmd::dcdx, w3adatmd::dcdy, w3adatmd::dcxdx, w3adatmd::dcxdy, w3adatmd::dcydx, w3adatmd::dcydy, w3adatmd::dddx, w3adatmd::dddy, w3gdatmd::dmin, w3gdatmd::dsip, w3gdatmd::dth, w3adatmd::dw, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::es2, w3gdatmd::esc, w3gdatmd::esin, w3gdatmd::fachfa, w3gdatmd::flck, w3gdatmd::flcth, w3idatmd::flcur, w3odatmd::iaproc, w3gdatmd::mapsf, w3gdatmd::mapwn, w3gdatmd::nk, w3gdatmd::nk2, w3gdatmd::nspec, w3gdatmd::nth, w3gdatmd::sig, w3servmd::strace(), and w3adatmd::wn.

Referenced by pdlib_w3profsmd::calcarray_jacobi_spectral_1(), and pdlib_w3profsmd::calcarray_jacobi_spectral_2().

◆ prop_refraction_pr3()

subroutine w3parall::prop_refraction_pr3 ( integer, intent(in)  IP,
integer, intent(in)  ISEA,
real, intent(in)  DTG,
real, dimension(nspec), intent(out)  CAD,
logical, intent(in)  DoLimiter 
)

Compute refraction part in matrix alternative approach.

Parameters
[in]IP
[in]ISEA
[in]DTG
[out]CAD
[in]DoLimiter
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 449 of file w3parall.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-June-2018 |
458  !/ +-----------------------------------+
459  !/
460  !/ 01-June-2018 : Origination. ( version 6.04 )
461  !/
462  ! 1. Purpose : Compute refraction part in matrix alternative approach
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  USE constants, ONLY : lpdlib
497  USE w3gdatmd, ONLY: nk, nk2, nth, nspec, sig, dsip, ecos, esin, &
498  ec2, esc, es2, fachfa, mapwn, flcth, flck, &
500  USE w3adatmd, ONLY: cg, wn, dcxdx, dcxdy, dcydx, dcydy, dddx, &
501  dddy, dw
502  USE w3idatmd, ONLY: flcur
503  USE w3odatmd, only : iaproc
504  IMPLICIT NONE
505  !/
506  !/ ------------------------------------------------------------------- /
507  !/ Parameter list
508  !/
509  !/ ------------------------------------------------------------------- /
510  !/ Local parameters
511  !/
512 #ifdef W3_S
513  INTEGER, SAVE :: IENT = 0
514 #endif
515  REAL, intent(out) :: CAD(NSPEC)
516  INTEGER, intent(in) :: ISEA, IP
517  REAL, intent(in) :: DTG
518  logical, intent(in) :: DoLimiter
519  INTEGER :: ISP, IK, ITH, IX, IY
520  REAL :: FRK(NK), FRG(NK), DSDD(0:NK+1)
521  REAL :: FACTH, DCXY, DCYX, DCXXYY, DTTST
522  REAL :: eDCXDX, eDCXDY, eDCYDX, eDCYDY, eDDDX, eDDDY, eCTHG0
523  REAL :: VCFLT(NSPEC), DEPTH, FDG, CG1(0:NK+1), WN1(0:NK+1)
524  REAL :: FDDMAX, CFLTHMAX, VELNOFILT, CTMAX_eff
525 #ifdef W3_S
526  CALL strace (ient, 'PROP_REFRACTION_PR3')
527 #endif
528 
529  ix = mapsf(isea,1)
530  iy = mapsf(isea,2)
531  edddx=dddx(1,ip)
532  edddy=dddy(1,ip)
533  ecthg0 = cthg0s(isea)
534  facth = 1.0 / dth
535  !
536  fdg = facth * ecthg0
537  depth = max( dmin , dw(isea) )
538  DO ik=0, nk+1
539  IF ( depth*wn(ik,isea) .LT. 5. ) THEN
540  dsdd(ik) = max( 0. , cg(ik,isea)*wn(ik,isea)-0.5*sig(ik) ) / depth
541  ELSE
542  dsdd(ik) = 0.
543  END IF
544  END DO
545  DO ik=1, nk
546  frk(ik) = facth * dsdd(ik) / wn(ik,isea)
547  frg(ik) = fdg * cg(ik,isea)
548  END DO
549  IF (flcur) THEN
550  edcxdx = dcxdx(1,ip)
551  edcxdy = dcxdy(1,ip)
552  edcydx = dcydx(1,ip)
553  edcydy = dcydy(1,ip)
554  dcyx = facth * edcydx
555  dcxxyy = facth * ( edcxdx - edcydy )
556  dcxy = facth * edcxdy
557  DO isp=1, nspec
558  vcflt(isp) = es2(isp)*dcyx + esc(isp)*dcxxyy - ec2(isp)*dcxy
559  END DO
560  ELSE
561  vcflt=0
562  END IF
563  !
564 #ifdef W3_REFRX
565  ! 3.c @C/@x refraction and great-circle propagation
566  DO ik=1, nk
567  frk(ik) = facth * cg(ik,isea) * wn(ik,isea) / sig(ik)
568  END DO
569 #endif
570  !
571  ctmax_eff=ctmax/dtg
572  DO isp=1, nspec
573  velnofilt = vcflt(isp) &
574  + frg(mapwn(isp)) * ecos(isp) &
575  + frk(mapwn(isp)) * (esin(isp)*edddx - ecos(isp)*edddy)
576  !
577  ! Puts filtering on total velocity (including currents and great circle effects)
578  ! the filtering limits VCFLT to be less than CTMAX
579  ! this modification was proposed by F. Ardhuin 2011/03/06
580  !
581  IF (dolimiter) THEN
582  vcflt(isp)=sign(min(abs(velnofilt),ctmax_eff),velnofilt)
583  ELSE
584  vcflt(isp)=velnofilt
585  END IF
586  END DO
587  DO isp=1,nspec
588  cad(isp)=dble(vcflt(isp))
589  END DO
590  !/
591  !/ End of JACOBI_INIT ------------------------------------------------ /
592  !/

References w3adatmd::cg, w3gdatmd::cthg0s, w3gdatmd::ctmax, w3adatmd::dcxdx, w3adatmd::dcxdy, w3adatmd::dcydx, w3adatmd::dcydy, w3adatmd::dddx, w3adatmd::dddy, w3gdatmd::dmin, w3gdatmd::dsip, w3gdatmd::dth, w3adatmd::dw, w3gdatmd::ec2, w3gdatmd::ecos, w3gdatmd::es2, w3gdatmd::esc, w3gdatmd::esin, w3gdatmd::fachfa, w3gdatmd::flck, w3gdatmd::flcth, w3idatmd::flcur, w3odatmd::iaproc, constants::lpdlib, w3gdatmd::mapsf, w3gdatmd::mapwn, w3gdatmd::nk, w3gdatmd::nk2, w3gdatmd::nspec, w3gdatmd::nth, w3gdatmd::sig, w3servmd::strace(), and w3adatmd::wn.

Referenced by pdlib_w3profsmd::calcarray_jacobi_spectral_1(), and pdlib_w3profsmd::calcarray_jacobi_spectral_2().

◆ set_up_nseal_nsealm()

subroutine w3parall::set_up_nseal_nsealm ( integer, intent(out)  NSEALout,
integer, intent(out)  NSEALMout 
)

Setup NSEAL, NSEALM in context of PDLIB.

Parameters
[out]NSEALout
[out]NSEALMout
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 1040 of file w3parall.F90.

1040  !/
1041  !/ +-----------------------------------+
1042  !/ | WAVEWATCH III NOAA/NCEP |
1043  !/ | |
1044  !/ | Aron Roland (BGS IT&E GmbH) |
1045  !/ | Mathieu Dutour-Sikiric (IRB) |
1046  !/ | |
1047  !/ | FORTRAN 90 |
1048  !/ | Last update : 01-June-2018 |
1049  !/ +-----------------------------------+
1050  !/
1051  !/ 01-June-2018 : Origination. ( version 6.04 )
1052  !/
1053  ! 1. Purpose : Setup nseal, nsealm in contect of pdlib
1054  ! 2. Method :
1055  ! 3. Parameters :
1056  !
1057  ! Parameter list
1058  ! ----------------------------------------------------------------
1059  ! ----------------------------------------------------------------
1060  !
1061  ! 4. Subroutines used :
1062  !
1063  ! Name Type Module Description
1064  ! ----------------------------------------------------------------
1065  ! STRACE Subr. W3SERVMD Subroutine tracing.
1066  ! ----------------------------------------------------------------
1067  !
1068  ! 5. Called by :
1069  !
1070  ! Name Type Module Description
1071  ! ----------------------------------------------------------------
1072  ! ----------------------------------------------------------------
1073  !
1074  ! 6. Error messages :
1075  ! 7. Remarks
1076  ! 8. Structure :
1077  ! 9. Switches :
1078  !
1079  ! !/S Enable subroutine tracing.
1080  !
1081  ! 10. Source code :
1082  !
1083  !/ ------------------------------------------------------------------- /
1084 #ifdef W3_S
1085  USE w3servmd, ONLY: strace
1086 #endif
1087  !/
1088  !/ ------------------------------------------------------------------- /
1089 #ifdef W3_PDLIB
1090  use yowdatapool, only: istatus
1091  use yownodepool, only: npa
1092  use yowrankmodule, only : rank
1093  USE w3gdatmd, ONLY: gtype, ungtype
1094 #endif
1095 #ifdef W3_MPI
1097 #endif
1098  USE constants, ONLY : lpdlib
1099  USE w3gdatmd, ONLY: nsea
1100  USE w3odatmd, ONLY: ntproc, naproc, iaproc
1101  IMPLICIT NONE
1102  INTEGER, intent(out) :: NSEALout, NSEALMout
1103  !/ Local parameters
1104  !/
1105 #ifdef W3_S
1106  INTEGER, SAVE :: IENT = 0
1107 #endif
1108  !/
1109  !/ ------------------------------------------------------------------- /
1110  !/
1111 #ifdef W3_S
1112  CALL strace (ient, 'SET_UP_NSEAL_NSEALM')
1113 #endif
1114 
1115 #ifdef W3_SHRD
1116  nsealout = nsea
1117  nsealmout = nsea
1118 #endif
1119  !
1120 #ifdef W3_DIST
1121  IF (.NOT. lpdlib ) THEN
1122  IF ( iaproc .LE. naproc ) THEN
1123  nsealout = 1 + (nsea-iaproc)/naproc
1124  ELSE
1125  nsealout = 0
1126  END IF
1127  nsealmout = 1 + (nsea-1)/naproc
1128  ELSE
1129 #endif
1130 #ifdef W3_PDLIB
1131  IF (gtype .eq. ungtype) THEN
1132  nsealout = pdlib_nseal
1133  nsealmout = pdlib_nsealm
1134  ELSE
1135  IF ( iaproc .LE. naproc ) THEN
1136  nsealout = 1 + (nsea-iaproc)/naproc
1137  ELSE
1138  nsealout = 0
1139  END IF
1140  nsealmout = 1 + (nsea-1)/naproc
1141  ENDIF
1142 #endif
1143 #ifdef W3_DIST
1144  ENDIF
1145 #endif
1146  !/
1147  !/ End of JACOBI_INIT ------------------------------------------------ /
1148  !/

References w3gdatmd::gtype, w3odatmd::iaproc, yowdatapool::istatus, constants::lpdlib, w3adatmd::mpi_comm_wave, w3adatmd::mpi_comm_wcmp, w3odatmd::naproc, yownodepool::npa, w3gdatmd::nsea, w3odatmd::ntproc, pdlib_nseal, pdlib_nsealm, yowrankmodule::rank, w3servmd::strace(), and w3gdatmd::ungtype.

Referenced by w3wdatmd::w3dimw(), and w3initmd::w3init().

◆ synchronize_global_array()

subroutine w3parall::synchronize_global_array ( double precision, dimension(nx), intent(inout)  TheVar)

Sync global array in context of PDLIB.

An array of size (NSEA) is send but only the (1:NSEAL) values are correct. The program synchonizes everything on all nodes.

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

Definition at line 1517 of file w3parall.F90.

1517  !/ ------------------------------------------------------------------- /
1518  !**********************************************************************
1519  !* An array of size (NSEA) is send but only the (1:NSEAL) values *
1520  !* are correct. The program synchonizes everything on all nodes. *
1521  !**********************************************************************
1522  !/
1523  !/ +-----------------------------------+
1524  !/ | WAVEWATCH III NOAA/NCEP |
1525  !/ | |
1526  !/ | Aron Roland (BGS IT&E GmbH) |
1527  !/ | Mathieu Dutour-Sikiric (IRB) |
1528  !/ | |
1529  !/ | FORTRAN 90 |
1530  !/ | Last update : 01-June-2018 |
1531  !/ +-----------------------------------+
1532  !/
1533  !/ 01-June-2018 : Origination. ( version 6.04 )
1534  !/
1535  ! 1. Purpose : Sync global array in context of pdlib
1536  ! 2. Method :
1537  ! An array of size (NSEA) is send but only the (1:NSEAL) values
1538  ! are correct. The program synchonizes everything on all nodes.
1539  ! 3. Parameters :
1540  !
1541  ! Parameter list
1542  ! ----------------------------------------------------------------
1543  ! ----------------------------------------------------------------
1544  !
1545  ! 4. Subroutines used :
1546  !
1547  ! Name Type Module Description
1548  ! ----------------------------------------------------------------
1549  ! STRACE Subr. W3SERVMD Subroutine tracing.
1550  ! ----------------------------------------------------------------
1551  !
1552  ! 5. Called by :
1553  !
1554  ! Name Type Module Description
1555  ! ----------------------------------------------------------------
1556  ! ----------------------------------------------------------------
1557  !
1558  ! 6. Error messages :
1559  ! 7. Remarks
1560  ! 8. Structure :
1561  ! 9. Switches :
1562  !
1563  ! !/S Enable subroutine tracing.
1564  !
1565  ! 10. Source code :
1566  !
1567  !/ ------------------------------------------------------------------- /
1568 #ifdef W3_S
1569  USE w3servmd, ONLY: strace
1570 #endif
1571  !
1572  USE w3gdatmd, ONLY: nseal, nsea, nx
1573 #ifdef W3_PDLIB
1574  USE w3odatmd, only : iaproc, naproc, ntproc
1575  USE w3adatmd, ONLY: mpi_comm_wcmp
1576  use yowdatapool, only: rtype, istatus
1577  USE yownodepool, only: npa
1578  use yownodepool, only: iplg
1579  use yowdatapool, only: rkind
1580 #endif
1581  IMPLICIT NONE
1582  !/ ------------------------------------------------------------------- /
1583  !/ Parameter list
1584  !/
1585  !/ ------------------------------------------------------------------- /
1586  !/ Local parameters
1587  !/
1588 #ifdef W3_S
1589  INTEGER, SAVE :: IENT = 0
1590 #endif
1591  !/
1592  !/ ------------------------------------------------------------------- /
1593  !/
1594 #ifdef W3_MPI
1595  include "mpif.h"
1596 #endif
1597  INTEGER ISEA, JSEA, Status(NX), rStatus(NX)
1598  INTEGER IPROC, I, ierr, IP, IX, IP_glob
1599 #ifdef W3_PDLIB
1600  REAL(rkind), intent(inout) :: TheVar(NX)
1601  REAL(rkind) :: rVect(NX)
1602 #else
1603  DOUBLE PRECISION, intent(inout) :: TheVar(NX)
1604  DOUBLE PRECISION :: rVect(NX)
1605 #endif
1606  status=0
1607 #ifdef W3_S
1608  CALL strace (ient, 'SYNCHRONIZE_GLOBAL_ARRAY')
1609 #endif
1610 #ifdef W3_PDLIB
1611  DO ip=1,npa
1612  ip_glob=iplg(ip)
1613  status(ip_glob)=1
1614  END DO
1615  IF (iaproc .eq. 1) THEN
1616  DO iproc=2,naproc
1617  CALL mpi_recv(rvect,nx,rtype, iproc-1, 19, mpi_comm_wcmp, istatus, ierr)
1618  CALL mpi_recv(rstatus,nx,mpi_integer, iproc-1, 23, mpi_comm_wcmp, istatus, ierr)
1619  DO i=1,nx
1620  IF (rstatus(i) .eq. 1) THEN
1621  thevar(i)=rvect(i)
1622  status(i)=1
1623  END IF
1624  END DO
1625  END DO
1626  DO iproc=2,naproc
1627  CALL mpi_send(thevar,nx,rtype, iproc-1, 29, mpi_comm_wcmp, ierr)
1628  END DO
1629  ELSE
1630  CALL mpi_send(thevar,nx,rtype, 0, 19, mpi_comm_wcmp, ierr)
1631  CALL mpi_send(status,nx,mpi_integer, 0, 23, mpi_comm_wcmp, ierr)
1632  CALL mpi_recv(thevar,nx,rtype, 0, 29, mpi_comm_wcmp, istatus, ierr)
1633  END IF
1634 #endif
1635  !/
1636  !/ End of JACOBI_INIT ------------------------------------------------ /
1637  !/

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

Referenced by w3wavset::trig_wave_setup_computation().

◆ synchronize_ipgl_etc_array()

subroutine w3parall::synchronize_ipgl_etc_array ( integer, intent(in)  IMOD,
logical, intent(in)  IsMulti 
)

Sync global local arrays.

Parameters
[in]IMOD
[in]IsMulti
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 916 of file w3parall.F90.

916  !/
917  !/ +-----------------------------------+
918  !/ | WAVEWATCH III NOAA/NCEP |
919  !/ | |
920  !/ | Aron Roland (BGS IT&E GmbH) |
921  !/ | Mathieu Dutour-Sikiric (IRB) |
922  !/ | |
923  !/ | FORTRAN 90 |
924  !/ | Last update : 01-June-2018 |
925  !/ +-----------------------------------+
926  !/
927  !/ 01-June-2018 : Origination. ( version 6.04 )
928  !/
929  ! 1. Purpose : Sync global local arrays
930  ! 2. Method :
931  ! All the process need to have IPGL_tot and IPGL_TO_PROC
932  ! This is especially the case for the output process.
933  ! So we need some painful exportation business
934  ! 3. Parameters :
935  !
936  ! Parameter list
937  ! ----------------------------------------------------------------
938  ! ----------------------------------------------------------------
939  !
940  ! 4. Subroutines used :
941  !
942  ! Name Type Module Description
943  ! ----------------------------------------------------------------
944  ! STRACE Subr. W3SERVMD Subroutine tracing.
945  ! ----------------------------------------------------------------
946  !
947  ! 5. Called by :
948  !
949  ! Name Type Module Description
950  ! ----------------------------------------------------------------
951  ! ----------------------------------------------------------------
952  !
953  ! 6. Error messages :
954  ! 7. Remarks
955  ! 8. Structure :
956  ! 9. Switches :
957  !
958  ! !/S Enable subroutine tracing.
959  !
960  ! 10. Source code :
961  !
962  !/ ------------------------------------------------------------------- /
963 #ifdef W3_S
964  USE w3servmd, ONLY: strace
965 #endif
966 #ifdef W3_PDLIB
967  USE yowdatapool, only: istatus
968  USE yownodepool, only: np_global
969  USE w3odatmd, ONLY: ntproc, naproc, iaproc
970  USE w3gdatmd, ONLY: mapsf, nsea
972  USE yowrankmodule, only : ipgl_to_proc, ipgl_tot
973  USE wmmdatmd, ONLY: mdatas
974 #endif
975  IMPLICIT NONE
976 #ifdef W3_PDLIB
977  include "mpif.h"
978 #endif
979  INTEGER, intent(in) :: IMOD
980  logical, intent(in) :: IsMulti
981 #ifdef W3_S
982  INTEGER, SAVE :: IENT = 0
983 #endif
984 #ifdef W3_PDLIB
985  INTEGER :: Iarr(1)
986  INTEGER :: ISEA, IP_glob
987  INTEGER :: IPROC, IERR_MPI, istat
988 #endif
989 
990 #ifdef W3_S
991  CALL strace (ient, 'SYNCHRONIZE_IPGL_ETC_ARRAY')
992 #endif
993 
994 #ifdef W3_PDLIB
995  IF (iaproc .le. naproc) THEN
996  IF (iaproc .eq. 1) THEN
997  iarr(1)=np_global
998  DO iproc=naproc+1,ntproc
999  CALL mpi_send(iarr,1,mpi_int, iproc-1, 37, mpi_comm_wave, ierr_mpi)
1000  END DO
1001  DO iproc=naproc+1,ntproc
1002  CALL mpi_send(ipgl_tot,np_global,mpi_int, iproc-1, 43, mpi_comm_wave, ierr_mpi)
1003  CALL mpi_send(ipgl_to_proc,np_global,mpi_int, iproc-1, 91, mpi_comm_wave, ierr_mpi)
1004  END DO
1005  END IF
1006  ELSE
1007  CALL mpi_recv(iarr,1,mpi_int, 0, 37, mpi_comm_wave, istatus, ierr_mpi)
1008  np_global=iarr(1)
1009  allocate(ipgl_tot(np_global), ipgl_to_proc(np_global), stat=istat)
1010  CALL mpi_recv(ipgl_tot,np_global,mpi_int, 0, 43, mpi_comm_wave, istatus, ierr_mpi)
1011  CALL mpi_recv(ipgl_to_proc,np_global,mpi_int, 0, 91, mpi_comm_wave, istatus, ierr_mpi)
1012  END IF
1013  IF (ismulti) THEN
1014  WRITE(*,*) ' Before allocation of MDATAS % SEA_IPGL, SEA_IPGL_TO_PROC : IMOD=', imod, ' NSEA=', nsea
1015  ALLOCATE(mdatas(imod)%SEA_IPGL(nsea), mdatas(imod)%SEA_IPGL_TO_PROC(nsea), stat=istat)
1016  !CHECK_ALLOC_STATUS ( ISTAT )
1017  DO isea=1,nsea
1018  ip_glob = mapsf(isea, 1)
1019  mdatas(imod)%SEA_IPGL(isea) = ipgl_tot(ip_glob)
1020  mdatas(imod)%SEA_IPGL_TO_PROC(isea) = ipgl_to_proc(ip_glob)
1021  END DO
1022  END IF
1023 #endif
1024  !/
1025  !/ End of JACOBI_INIT ------------------------------------------------ /
1026  !/

References w3odatmd::iaproc, include(), yowrankmodule::ipgl_to_proc, yowrankmodule::ipgl_tot, yowdatapool::istatus, w3gdatmd::mapsf, wmmdatmd::mdatas, w3adatmd::mpi_comm_wave, w3adatmd::mpi_comm_wcmp, w3odatmd::naproc, yownodepool::np_global, w3gdatmd::nsea, w3odatmd::ntproc, and w3servmd::strace().

Referenced by w3initmd::w3init().

◆ wav_my_wtime()

subroutine w3parall::wav_my_wtime ( real(8), intent(out)  eTime)

NA.

Parameters
[out]eTime
Author
Aron Roland
Mathieu Dutour-Sikiric
Date
01-Jun-2018

Definition at line 110 of file w3parall.F90.

110  !/ ------------------------------------------------------------------- /
111  !/
112  !/ +-----------------------------------+
113  !/ | WAVEWATCH III NOAA/NCEP |
114  !/ | |
115  !/ | Aron Roland (BGS IT&E GmbH) |
116  !/ | Mathieu Dutour-Sikiric (IRB) |
117  !/ | |
118  !/ | FORTRAN 90 |
119  !/ | Last update : 01-June-2018 |
120  !/ +-----------------------------------+
121  !/
122  !/ 01-June-2018 : Origination. ( version 6.04 )
123  !/
124  ! 1. Purpose :
125  ! 2. Method :
126  ! 3. Parameters :
127  !
128  ! Parameter list
129  ! ----------------------------------------------------------------
130  ! ----------------------------------------------------------------
131  !
132  ! 4. Subroutines used :
133  !
134  ! Name Type Module Description
135  ! ----------------------------------------------------------------
136  ! STRACE Subr. W3SERVMD Subroutine tracing.
137  ! ----------------------------------------------------------------
138  !
139  ! 5. Called by :
140  !
141  ! Name Type Module Description
142  ! ----------------------------------------------------------------
143  ! ----------------------------------------------------------------
144  !
145  ! 6. Error messages :
146  ! 7. Remarks
147  ! 8. Structure :
148  ! 9. Switches :
149  !
150  ! !/S Enable subroutine tracing.
151  !
152  ! 10. Source code :
153  !
154  !/ ------------------------------------------------------------------- /
155 #ifdef W3_S
156  USE w3servmd, ONLY: strace
157 #endif
158  !/
159  !/ ------------------------------------------------------------------- /
160  !/ Parameter list
161  !/
162  !/ ------------------------------------------------------------------- /
163  !/ Local parameters
164  !/
165  IMPLICIT NONE
166 #ifdef W3_S
167  INTEGER, SAVE :: IENT = 0
168 #endif
169  INTEGER mpimode
170  REAL(8), intent(out) :: eTime
171 #ifdef W3_MPI
172  REAL(8) mpi_wtime
173 #endif
174  mpimode=0
175 #ifdef W3_MPI
176  mpimode=1
177  etime=mpi_wtime()
178 #endif
179 #ifdef W3_S
180  CALL strace (ient, 'WAV_MY_WTIME')
181 #endif
182  IF (mpimode .eq. 0) THEN
183  CALL cpu_time(etime)
184  END IF
185  !/
186  !/ End of JACOBI_INIT ------------------------------------------------ /
187  !/

References w3servmd::strace().

Referenced by print_my_time().

Variable Documentation

◆ ient

integer, save w3parall::ient = 0

Definition at line 78 of file w3parall.F90.

78  INTEGER, SAVE :: IENT = 0

◆ imem

◆ isea_to_jsea

integer, dimension(:), allocatable w3parall::isea_to_jsea

◆ jx_to_jsea

integer, dimension(:), allocatable w3parall::jx_to_jsea

Definition at line 83 of file w3parall.F90.

83  INTEGER, ALLOCATABLE :: JX_TO_JSEA(:), ISEA_TO_JSEA(:)

Referenced by get_jsea_ibelong(), pdlib_w3profsmd::pdlib_init(), and pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block().

◆ listispnextdir

integer, dimension(:), allocatable w3parall::listispnextdir

Definition at line 86 of file w3parall.F90.

86  INTEGER, ALLOCATABLE :: ListISPnextDir(:), ListISPprevDir(:)

Referenced by pdlib_w3profsmd::block_solver_init(), and pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block().

◆ listispnextfreq

integer, dimension(:), allocatable w3parall::listispnextfreq

Definition at line 87 of file w3parall.F90.

87  INTEGER, ALLOCATABLE :: ListISPnextFreq(:), ListISPprevFreq(:)

Referenced by pdlib_w3profsmd::block_solver_init().

◆ listispprevdir

integer, dimension(:), allocatable w3parall::listispprevdir

◆ listispprevfreq

integer, dimension(:), allocatable w3parall::listispprevfreq

Definition at line 87 of file w3parall.F90.

Referenced by pdlib_w3profsmd::block_solver_init().

◆ lsloc

logical, parameter w3parall::lsloc = .true.

Definition at line 89 of file w3parall.F90.

89  LOGICAL, PARAMETER :: LSLOC = .true.

Referenced by pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block(), w3wdatmd::w3dimw(), w3srcemd::w3srce(), and w3wavemd::w3wave().

◆ onesixth

◆ onethird

real, parameter w3parall::onethird = 1.0d0/3.0d0

Definition at line 93 of file w3parall.F90.

93  REAL, PARAMETER :: ONETHIRD = 1.0d0/3.0d0

Referenced by pdlib_w3profsmd::calcarray_jacobi3(), and pdlib_w3profsmd::calcarray_jacobi4().

◆ pdlib_nseal

integer w3parall::pdlib_nseal

Definition at line 82 of file w3parall.F90.

82  INTEGER :: PDLIB_NSEAL, PDLIB_NSEALM

Referenced by pdlib_w3profsmd::pdlib_init(), set_up_nseal_nsealm(), and w3wavemd::w3wave().

◆ pdlib_nsealm

integer w3parall::pdlib_nsealm

Definition at line 82 of file w3parall.F90.

Referenced by pdlib_w3profsmd::pdlib_init(), set_up_nseal_nsealm(), and w3wavemd::w3wave().

◆ thr

◆ thr8

real*8, parameter w3parall::thr8 = TINY(1.d0)

Definition at line 96 of file w3parall.F90.

96  real*8, PARAMETER :: thr8 = tiny(1.d0)

Referenced by pdlib_w3profsmd::pdlib_jacobi_gauss_seidel_block().

◆ zero

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3gdatmd::esc
real, dimension(:), pointer esc
Definition: w3gdatmd.F90:1234
w3gdatmd::nseal
integer, pointer nseal
Definition: w3gdatmd.F90:1097
include
cmake src_list cmake include(${CMAKE_CURRENT_SOURCE_DIR}/cmake/check_switches.cmake) check_switches("$
Definition: CMakeLists.txt:15
yowrankmodule::ipgl_tot
integer, dimension(:), allocatable, public ipgl_tot
Definition: yowrankModule.F90:69
w3gdatmd::dth
real, pointer dth
Definition: w3gdatmd.F90:1232
w3adatmd
Define data structures to set up wave model auxiliary data for several models simultaneously.
Definition: w3adatmd.F90:26
w3gdatmd::nspec
integer, pointer nspec
Definition: w3gdatmd.F90:1230
yowrankmodule::rank
type(t_rank), dimension(:), allocatable, public rank
Provides access to some information of all threads e.g.
Definition: yowrankModule.F90:68
w3adatmd::dcdx
real, dimension(:,:,:), pointer dcdx
Definition: w3adatmd.F90:629
w3gdatmd::ungtype
integer, parameter ungtype
Definition: w3gdatmd.F90:626
yowrankmodule::ipgl_to_proc
integer, dimension(:), allocatable, public ipgl_to_proc
Definition: yowrankModule.F90:69
w3gdatmd::dmin
real, pointer dmin
Definition: w3gdatmd.F90:1183
w3adatmd::cg
real, dimension(:,:), pointer cg
Definition: w3adatmd.F90:575
w3adatmd::dw
real, dimension(:), pointer dw
Definition: w3adatmd.F90:584
w3odatmd::ntproc
integer, pointer ntproc
Definition: w3odatmd.F90:457
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3gdatmd::flck
logical, pointer flck
Definition: w3gdatmd.F90:1217
w3odatmd::iaproc
integer, pointer iaproc
Definition: w3odatmd.F90:457
w3gdatmd::ecos
real, dimension(:), pointer ecos
Definition: w3gdatmd.F90:1234
yownodepool::iplg
integer, dimension(:), allocatable, public iplg
Node local to global mapping.
Definition: yownodepool.F90:116
yownodepool::npa
integer, public npa
number of ghost + resident nodes this partition holds
Definition: yownodepool.F90:99
w3gdatmd::dsip
real, dimension(:), pointer dsip
Definition: w3gdatmd.F90:1234
w3idatmd::flcur
logical, pointer flcur
Definition: w3idatmd.F90:261
yowdatapool::rkind
integer, parameter rkind
double precision.
Definition: yowdatapool.F90:47
w3adatmd::dcydx
real, dimension(:,:), pointer dcydx
Definition: w3adatmd.F90:627
scrip_constants::tiny
real(kind=scrip_r8), parameter tiny
Definition: scrip_constants.f:50
yownodepool::np_global
integer, public np_global
number of nodes, global
Definition: yownodepool.F90:89
w3gdatmd::fachfa
real, pointer fachfa
Definition: w3gdatmd.F90:1232
yowrankmodule::ipgl_npa
integer, dimension(:), allocatable, public ipgl_npa
Definition: yowrankModule.F90:70
w3gdatmd::es2
real, dimension(:), pointer es2
Definition: w3gdatmd.F90:1234
w3gdatmd::esin
real, dimension(:), pointer esin
Definition: w3gdatmd.F90:1234
scrip_constants::zero
real(kind=scrip_r8), parameter zero
Definition: scrip_constants.f:50
yownodepool::ipgl
integer, dimension(:), allocatable, public ipgl
Node global to local mapping np_global long.
Definition: yownodepool.F90:120
w3adatmd::dcdy
real, dimension(:,:,:), pointer dcdy
Definition: w3adatmd.F90:629
wmmdatmd::mdatas
type(mdata), dimension(:), allocatable, target mdatas
MDATAS.
Definition: wmmdatmd.F90:528
yownodepool
Has data that belong to nodes.
Definition: yownodepool.F90:39
constants::lpdlib
logical lpdlib
LPDLIB Logical for using the PDLIB library.
Definition: constants.F90:101
w3gdatmd::nk2
integer, pointer nk2
Definition: w3gdatmd.F90:1230
w3gdatmd::nsea
integer, pointer nsea
Definition: w3gdatmd.F90:1097
w3servmd
Definition: w3servmd.F90:3
w3adatmd::dcxdx
real, dimension(:,:), pointer dcxdx
Definition: w3adatmd.F90:627
yowrankmodule
Provides access to some information of all threads e.g.
Definition: yowrankModule.F90:44
yowdatapool::rtype
integer, save rtype
Definition: yowdatapool.F90:76
w3gdatmd::flcth
logical, pointer flcth
Definition: w3gdatmd.F90:1217
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3odatmd
Definition: w3odatmd.F90:3
w3adatmd::dddy
real, dimension(:,:), pointer dddy
Definition: w3adatmd.F90:627
w3adatmd::cy
real, dimension(:), pointer cy
Definition: w3adatmd.F90:584
w3gdatmd::mapsf
integer, dimension(:,:), pointer mapsf
Definition: w3gdatmd.F90:1163
w3odatmd::naproc
integer, pointer naproc
Definition: w3odatmd.F90:457
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
w3gdatmd::cthg0s
real, dimension(:), pointer cthg0s
Definition: w3gdatmd.F90:1198
w3adatmd::wn
real, dimension(:,:), pointer wn
Definition: w3adatmd.F90:575
w3adatmd::dcydy
real, dimension(:,:), pointer dcydy
Definition: w3adatmd.F90:627
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3gdatmd::gtype
integer, pointer gtype
Definition: w3gdatmd.F90:1094
w3idatmd
Define data structures to set up wave model input data for several models simultaneously.
Definition: w3idatmd.F90:16
w3gdatmd::mapwn
integer, dimension(:), pointer mapwn
Definition: w3gdatmd.F90:1231
w3adatmd::mpi_comm_wave
integer, pointer mpi_comm_wave
Definition: w3adatmd.F90:676
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
wmmdatmd
Define data structures to set up wave model dynamic data for several models simultaneously.
Definition: wmmdatmd.F90:16
w3adatmd::dcxdy
real, dimension(:,:), pointer dcxdy
Definition: w3adatmd.F90:627
w3odatmd::outpts
type(output), dimension(:), allocatable, target outpts
Definition: w3odatmd.F90:452
w3adatmd::dddx
real, dimension(:,:), pointer dddx
Definition: w3adatmd.F90:627
w3adatmd::cx
real, dimension(:), pointer cx
Definition: w3adatmd.F90:584
w3gdatmd::nx
integer, pointer nx
Definition: w3gdatmd.F90:1097
w3gdatmd::ctmax
real, pointer ctmax
Definition: w3gdatmd.F90:1183
w3gdatmd::ec2
real, dimension(:), pointer ec2
Definition: w3gdatmd.F90:1234
w3adatmd::mpi_comm_wcmp
integer, pointer mpi_comm_wcmp
Definition: w3adatmd.F90:676