WAVEWATCH III  beta 0.0.1
w3snl1md Module Reference

Bundles routines to calculate nonlinear wave-wave interactions according to the Discrete Interaction Approximation (DIA) of Hasselmann et al. More...

Functions/Subroutines

subroutine w3snl1 (A, CG, KDMEAN, S, D)
 Calculate nonlinear interactions and the diagonal term of its derivative. More...
 
subroutine insnl1 (IMOD)
 Preprocessing for nonlinear interactions (weights). More...
 
subroutine w3snlgqm (A, CG, WN, DEPTH, TSTOTn, TSDERn)
 
double precision function couple (XK1, YK1, XK2, YK2, XK3, YK3, XK4, YK4)
 
subroutine gauleg (W_LEG, X_LEG, NPOIN)
 
subroutine f1f1f1 (F1SF, NF1, IQ_OM1)
 
subroutine insnlgqm
 

Variables

integer nconf
 
integer, dimension(:,:,:), allocatable k_if2
 
integer, dimension(:,:,:), allocatable k_if3
 
integer, dimension(:,:,:), allocatable k_1p2p
 
integer, dimension(:,:,:), allocatable k_1p3m
 
integer, dimension(:,:,:), allocatable k_1p2m
 
integer, dimension(:,:,:), allocatable k_1p3p
 
integer, dimension(:,:,:), allocatable k_1m2p
 
integer, dimension(:,:,:), allocatable k_1m3m
 
integer, dimension(:,:,:), allocatable k_1m2m
 
integer, dimension(:,:,:), allocatable k_1m3p
 
integer, dimension(:), allocatable f_poin
 
integer, dimension(:), allocatable t_poin
 
integer, dimension(:), allocatable k_if1
 
integer, dimension(:,:), allocatable k_1p
 
integer, dimension(:,:), allocatable k_1m
 
integer, dimension(:,:), allocatable idconf
 
double precision, dimension(:), allocatable f_coef
 
double precision, dimension(:), allocatable f_proj
 
double precision, dimension(:), allocatable tb_sca
 
double precision, dimension(:), allocatable tb_v14
 
double precision, dimension(:,:,:), allocatable tb_v24
 
double precision, dimension(:,:,:), allocatable tb_v34
 
double precision, dimension(:,:,:), allocatable tb_tpm
 
double precision, dimension(:,:,:), allocatable tb_tmp
 
double precision, dimension(:,:,:), allocatable tb_fac
 

Detailed Description

Bundles routines to calculate nonlinear wave-wave interactions according to the Discrete Interaction Approximation (DIA) of Hasselmann et al.

(JPO, 1985).

Author
H. L. Tolman
Date
03-Sep-2012

Function/Subroutine Documentation

◆ couple()

double precision function w3snl1md::couple ( double precision, intent(in)  XK1,
double precision, intent(in)  YK1,
double precision, intent(in)  XK2,
double precision, intent(in)  YK2,
double precision, intent(in)  XK3,
double precision, intent(in)  YK3,
double precision, intent(in)  XK4,
double precision, intent(in)  YK4 
)

Definition at line 1184 of file w3snl1md.F90.

1184  !/
1185  !/ +-----------------------------------+
1186  !/ | WAVEWATCH III NOAA/NCEP |
1187  !/ | M. Benoit & E. Gagnaire-Renou |
1188  !/ | Last update : 20-Nov-2022 |
1189  !/ +-----------------------------------+
1190  !/
1191  !/ 19-Nov-2022 : Transfer from TOMAWAC code ( version 7.xx )
1192  !/
1193  ! 1. Purpose :
1194  !
1195  ! Computes the 4-wave coupling coefficient used in Snl4
1196  !
1197  ! 2. Method :
1198  !
1199  ! Uses theoretical expression by Webb (1978)
1200  !
1201  ! 3. Parameters :
1202  !
1203  ! Parameter list
1204  ! ----------------------------------------------------------------
1205  ! XK1 Real I x component of k1 wavenumber ...
1206  ! ----------------------------------------------------------------
1207  !
1208  ! 5. Called by :
1209  !
1210  ! Name Type Module Description
1211  ! ----------------------------------------------------------------
1212  ! INNSLGQM Subr. W3SNL2 Prepares source term integration.
1213  ! ----------------------------------------------------------------
1214  !
1215  ! 6. Error messages :
1216  !
1217  ! None.
1218  !
1219  ! 10. Source code :
1220  !
1221  !/ ------------------------------------------------------------------- /
1222  USE constants, ONLY: grav
1223  !
1224  IMPLICIT NONE
1225 
1226  DOUBLE PRECISION, INTENT(IN) :: XK1 , YK1 , XK2 , YK2
1227  DOUBLE PRECISION, INTENT(IN) :: XK3 , YK3
1228  DOUBLE PRECISION, INTENT(IN) :: XK4 , YK4
1229  DOUBLE PRECISION COUPLE
1230  !
1231  !.....LOCAL VARIABLES
1232  ! """"""""""""""""""
1233  DOUBLE PRECISION RK1 , RK2 , RK3 , RK4 , WK1 , WK2
1234  DOUBLE PRECISION WK3 , WK4 , S12 , S13 , S14 , S23
1235  DOUBLE PRECISION S24 , S34 , W1P2 , Q12 , W1M3 , Q13
1236  DOUBLE PRECISION W1M4 , Q14 , DDD , COEF , DENO13, NUME13
1237  DOUBLE PRECISION DENO14, NUME14, ZERO, PI
1238 
1239  !
1240  pi = acos(-1.)
1241  coef=pi*grav*grav/4.d0
1242  zero=1.d-10
1243  !
1244  rk1=sqrt(xk1*xk1+yk1*yk1)
1245  rk2=sqrt(xk2*xk2+yk2*yk2)
1246  rk3=sqrt(xk3*xk3+yk3*yk3)
1247  rk4=sqrt(xk4*xk4+yk4*yk4)
1248  !
1249  wk1=sqrt(rk1)
1250  wk2=sqrt(rk2)
1251  wk3=sqrt(rk3)
1252  wk4=sqrt(rk4)
1253  !
1254  s12=xk1*xk2+yk1*yk2
1255  s13=xk1*xk3+yk1*yk3
1256  s14=xk1*xk4+yk1*yk4
1257  s23=xk2*xk3+yk2*yk3
1258  s24=xk2*xk4+yk2*yk4
1259  s34=xk3*xk4+yk3*yk4
1260  !
1261  w1p2=sqrt((xk1+xk2)*(xk1+xk2)+(yk1+yk2)*(yk1+yk2))
1262  w1m3=sqrt((xk1-xk3)*(xk1-xk3)+(yk1-yk3)*(yk1-yk3))
1263  w1m4=sqrt((xk1-xk4)*(xk1-xk4)+(yk1-yk4)*(yk1-yk4))
1264  q12=(wk1+wk2)*(wk1+wk2)
1265  q13=(wk1-wk3)*(wk1-wk3)
1266  q14=(wk1-wk4)*(wk1-wk4)
1267  !
1268  !.....COMPUTES THE D COEFFICIENT OF WEBB (1978)
1269  ! """"""""""""""""""""""""""""""""""""""
1270  ddd=2.00d0*q12*(rk1*rk2-s12)*(rk3*rk4-s34)/(w1p2-q12) &
1271  +0.50d0*(s12*s34+s13*s24+s14*s23) &
1272  +0.25d0*(s13+s24)*q13*q13 &
1273  -0.25d0*(s12+s34)*q12*q12 &
1274  +0.25d0*(s14+s23)*q14*q14 &
1275  +2.50d0*rk1*rk2*rk3*rk4 &
1276  +q12*q13*q14*(rk1+rk2+rk3+rk4)
1277 
1278  deno13=w1m3-q13
1279  nume13=2.00d0*q13*(rk1*rk3+s13)*(rk2*rk4+s24)
1280  IF (abs(deno13).LT.zero) THEN
1281  IF (abs(nume13).LT.zero) THEN
1282  WRITE(*,*) 'W3SNL2 error for coupling coefficient : (1-3) 0/0 !'
1283  ELSE
1284  WRITE(*,*) 'W3SNL2 error for coupling coefficient : (1-3) inifinte value'
1285  ENDIF
1286  WRITE(*,*) 'W3SNL2 error for coupling coefficient : (1-3) term not used'
1287  ELSE
1288  ddd=ddd+nume13/deno13
1289  ENDIF
1290  deno14=w1m4-q14
1291  nume14=2.00d0*q14*(rk1*rk4+s14)*(rk2*rk3+s23)
1292  IF (abs(deno14).LT.zero) THEN
1293  IF (abs(nume14).LT.zero) THEN
1294  WRITE(*,*) 'W3SNL2 error for coupling coefficient : (1-4) 0/0 !'
1295  ELSE
1296  WRITE(*,*) 'W3SNL2 error for coupling coefficient : (1-4) inifinte value'
1297  ENDIF
1298  WRITE(*,*) 'W3SNL2 error for coupling coefficient : (1-4) term not used'
1299  ELSE
1300  ddd=ddd+nume14/deno14
1301  ENDIF
1302 
1303  couple=coef*ddd*ddd/(wk1*wk2*wk3*wk4)
1304  ! RETURN

References constants::grav.

Referenced by insnlgqm().

◆ f1f1f1()

subroutine w3snl1md::f1f1f1 ( double precision, dimension(*), intent(inout)  F1SF,
integer, intent(inout)  NF1,
integer, intent(in)  IQ_OM1 
)

Definition at line 1347 of file w3snl1md.F90.

1347  ! TOMAWAC V6P3 15/06/2011
1348  !***********************************************************************
1349  !
1350  !brief SUBROUTINE CALLED BY PRENL3
1351  !+ COMPUTES VALUES OF RATIO F1/F AS FUNCTION OF THE IQ_OM1
1352  !+ INDICATOR
1353  !
1354  !history E. GAGNAIRE-RENOU
1355  !+ 04/2011
1356  !+ V6P1
1357  !+ CREATED
1358  !
1359  !history G.MATTAROLO (EDF - LNHE)
1360  !+ 15/06/2011
1361  !+ V6P1
1362  !+ Translation of French names of the variables in argument
1363  !
1364  !history E. GAGNAIRE-RENOU
1365  !+ 12/03/2013
1366  !+ V6P3
1367  !+ Better formatted: WRITE(LU,*), etc.
1368  !/ ------------------------------------------------------------------- /
1369  IMPLICIT NONE
1370  INTEGER, INTENT(IN) :: IQ_OM1
1371  INTEGER, INTENT(INOUT) :: NF1
1372  DOUBLE PRECISION, INTENT(INOUT) :: F1SF(*)
1373  !
1374  INTEGER I,M
1375  DOUBLE PRECISION RAISON
1376  !
1377  IF(iq_om1.EQ.1) THEN
1378  IF(nf1.NE.14) THEN
1379  WRITE(*,*) '#1 Incorrect value for NF1',nf1
1380  ENDIF
1381  f1sf( 1)=0.30d0
1382  f1sf( 2)=0.40d0
1383  f1sf( 3)=0.50d0
1384  f1sf( 4)=0.60d0
1385  f1sf( 5)=0.70d0
1386  f1sf( 6)=0.80d0
1387  f1sf( 7)=0.90d0
1388  f1sf( 8)=1.00d0
1389  f1sf( 9)=1.11d0
1390  f1sf(10)=1.25d0
1391  f1sf(11)=1.42d0
1392  f1sf(12)=1.67d0
1393  f1sf(13)=2.00d0
1394  f1sf(14)=2.50d0
1395  f1sf(15)=3.30d0
1396  ELSEIF(iq_om1.EQ.2) THEN
1397  IF (nf1.NE.26) THEN
1398  WRITE(*,*) '#2 Incorrect value for NF1', nf1
1399  ENDIF
1400  f1sf( 1)=0.32d0
1401  f1sf( 2)=0.35d0
1402  f1sf( 3)=0.39d0
1403  f1sf( 4)=0.44d0
1404  f1sf( 5)=0.50d0
1405  f1sf( 6)=0.56d0
1406  f1sf( 7)=0.63d0
1407  f1sf( 8)=0.70d0
1408  f1sf( 9)=0.78d0
1409  f1sf(10)=0.86d0
1410  f1sf(11)=0.92d0
1411  f1sf(12)=0.97d0
1412  f1sf(13)=1.00d0
1413  f1sf(14)=1.03d0
1414  f1sf(15)=1.08d0
1415  f1sf(16)=1.13d0
1416  f1sf(17)=1.20d0
1417  f1sf(18)=1.28d0
1418  f1sf(19)=1.37d0
1419  f1sf(20)=1.48d0
1420  f1sf(21)=1.50d0
1421  f1sf(22)=1.65d0
1422  f1sf(23)=1.85d0
1423  f1sf(24)=2.10d0
1424  f1sf(25)=2.40d0
1425  f1sf(26)=2.70d0
1426  f1sf(27)=3.20d0
1427  ELSEIF(iq_om1.EQ.3) THEN
1428  IF(nf1.NE.11) THEN
1429  WRITE(*,*) 'Incorrect value for NF1', nf1
1430  ENDIF
1431  f1sf( 1)=0.30d0
1432  f1sf( 2)=0.48d0
1433  f1sf( 3)=0.64d0
1434  f1sf( 4)=0.78d0
1435  f1sf( 5)=0.90d0
1436  f1sf( 6)=1.00d0
1437  f1sf( 7)=1.12d0
1438  f1sf( 8)=1.28d0
1439  f1sf( 9)=1.50d0
1440  f1sf(10)=1.80d0
1441  f1sf(11)=2.40d0
1442  f1sf(12)=3.40d0
1443  ELSEIF(iq_om1.EQ.4) THEN
1444  IF(nf1.NE.40) THEN
1445  WRITE(*,*) 'Incorrect value for NF1', nf1
1446  ENDIF
1447  nf1=20
1448  m=10
1449  raison=9.d0**(1.d0/dble(nf1))
1450  f1sf(m+1)=1.0d0/3.0d0
1451  nf1=2*m+nf1
1452  DO i=m+2,nf1+1
1453  f1sf(i)=f1sf(i-1)*raison
1454  ENDDO
1455  DO i=m,1,-1
1456  f1sf(i)=f1sf(i+1)/raison
1457  ENDDO
1458  ELSEIF(iq_om1.EQ.5) THEN
1459  raison=9.d0**(1.d0/dble(nf1))
1460  f1sf(1)=1.d0/3.d0
1461  DO i=2,nf1+1
1462  f1sf(i)=f1sf(i-1)*raison
1463  ENDDO
1464  ELSEIF(iq_om1.EQ.6) THEN
1465  raison=(3.d0-1.d0/3.d0)/dble(nf1)
1466  f1sf(1)=1.d0/3.d0
1467  DO i=2,nf1+1
1468  f1sf(i)=f1sf(i-1)+raison
1469  ENDDO
1470  ELSEIF(iq_om1.EQ.7) THEN
1471  IF(nf1.NE.20) THEN
1472  WRITE(*,*) 'Incorrect value for NF1', nf1
1473  ENDIF
1474  f1sf( 1)=1.d0/3.d0
1475  f1sf( 2)=0.40d0
1476  f1sf( 3)=0.46d0
1477  f1sf( 4)=0.52d0
1478  f1sf( 5)=0.60d0
1479  f1sf( 6)=0.70d0
1480  f1sf( 7)=0.79d0
1481  f1sf( 8)=0.86d0
1482  f1sf( 9)=0.92d0
1483  f1sf(10)=0.97d0
1484  f1sf(11)=1.00d0
1485  f1sf(12)=1.04d0
1486  f1sf(13)=1.10d0
1487  f1sf(14)=1.18d0
1488  f1sf(15)=1.28d0
1489  f1sf(16)=1.42d0
1490  f1sf(17)=1.60d0
1491  f1sf(18)=1.84d0
1492  f1sf(19)=2.14d0
1493  f1sf(20)=2.52d0
1494  f1sf(21)=3.00d0
1495  ENDIF
1496  !

Referenced by insnlgqm().

◆ gauleg()

subroutine w3snl1md::gauleg ( double precision, dimension(npoin), intent(inout)  W_LEG,
double precision, dimension(npoin), intent(inout)  X_LEG,
integer, intent(in)  NPOIN 
)

Definition at line 1309 of file w3snl1md.F90.

1309  !/ ------------------------------------------------------------------- /
1310  !.....VARIABLES IN ARGUMENT
1311  ! """"""""""""""""""""
1312  IMPLICIT NONE
1313  INTEGER , INTENT(IN) :: NPOIN
1314  DOUBLE PRECISION ,INTENT(INOUT) :: W_LEG(NPOIN) , X_LEG(NPOIN)
1315  !
1316  !.....LOCAL VARIABLES
1317  ! """""""""""""""""
1318  INTEGER I, M, J
1319  DOUBLE PRECISION EPS, Z, P1, P2, P3, PP, Z1, PI
1320  parameter(eps=3.d-14)
1321  !
1322  pi = acos(-1.)
1323  m=(npoin+1)/2
1324  DO i=1,m
1325  z=cos(pi*(dble(i)-0.25d0)/(dble(npoin)+0.5d0))
1326 1 CONTINUE
1327  p1=1.0d0
1328  p2=0.0d0
1329  DO j=1,npoin
1330  p3=p2
1331  p2=p1
1332  p1=((2.d0*dble(j)-1.d0)*z*p2-(dble(j)-1.d0)*p3)/dble(j)
1333  ENDDO
1334  pp=dble(npoin)*(z*p1-p2)/(z*z-1.d0)
1335  z1=z
1336  z=z-p1/pp
1337  IF (abs(z-z1).GT.eps) GOTO 1
1338  x_leg(i)=-z
1339  x_leg(npoin+1-i)=z
1340  w_leg(i)=2.d0/((1.d0-z**2)*pp**2)
1341  w_leg(npoin+1-i)=w_leg(i)
1342  ENDDO

Referenced by insnlgqm().

◆ insnl1()

subroutine w3snl1md::insnl1 ( integer, intent(in)  IMOD)

Preprocessing for nonlinear interactions (weights).

Parameters
[in]IMODModel number.
Author
H. L. Tolman
Date
24-Dec-2004

Definition at line 483 of file w3snl1md.F90.

483  !/
484  !/ +-----------------------------------+
485  !/ | WAVEWATCH III NOAA/NCEP |
486  !/ | H. L. Tolman |
487  !/ | FORTRAN 90 |
488  !/ | Last update : 24-Dec-2004 |
489  !/ +-----------------------------------+
490  !/
491  !/ 19-Oct-1998 : Final FORTRAN 77 ( version 1.18 )
492  !/ 04-Feb-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
493  !/ 09-May-2002 : Switch clean up. ( version 2.21 )
494  !/ 24-Dec-2004 : Multiple grid version. ( version 3.06 )
495  !/
496  ! 1. Purpose :
497  !
498  ! Preprocessing for nonlinear interactions (weights).
499  !
500  ! 2. Method :
501  !
502  ! See W3SNL1.
503  !
504  ! 3. Parameters :
505  !
506  ! Parameter list
507  ! ----------------------------------------------------------------
508  ! IMOD Int. I Model number.
509  ! ----------------------------------------------------------------
510  !
511  ! Local variables
512  ! ----------------------------------------------------------------
513  ! ITHxn Real Directional indices. (relative)
514  ! IFRxn Real Frequency indices. (relative)
515  ! IT1 R.A. Directional indices. (1-D)
516  ! IFn R.A. Frequency indices. (1-D)
517  ! ----------------------------------------------------------------
518  !
519  ! 4. Subroutines used :
520  !
521  ! Name Type Module Description
522  ! ----------------------------------------------------------------
523  ! STRACE Subr. W3SERVMD Subroutine tracing.
524  ! ----------------------------------------------------------------
525  !
526  ! 5. Called by :
527  !
528  ! Name Type Module Description
529  ! ----------------------------------------------------------------
530  ! W3IOGR Subr. W3IOGRMD Model definition file processing.
531  ! ----------------------------------------------------------------
532  !
533  ! 6. Error messages :
534  !
535  ! - Check on array dimensions for local arrays in W3SNL.
536  !
537  ! 7. Remarks :
538  !
539  ! - Test output is generated through W3IOGR.
540  ! - No testing of IMOD ir resetting of pointers.
541  !
542  ! 8. Structure :
543  !
544  ! - See source code.
545  !
546  ! 9. Switches :
547  !
548  ! !/S Enable subroutine tracing.
549  !
550  ! 10. Source code :
551  !
552  !/ ------------------------------------------------------------------- /
553  USE constants
554  USE w3gdatmd, ONLY: nk, nth, nspec, dth, xfr, sig, lam
555  USE w3adatmd, ONLY: w3dmnl
556  USE w3adatmd, ONLY: nfr, nfrhgh, nfrchg, nspecx, nspecy, &
557  ip11, ip12, ip13, ip14, im11, im12, im13, im14, &
558  ip21, ip22, ip23, ip24, im21, im22, im23, im24, &
559  ic11, ic12, ic21, ic22, ic31, ic32, ic41, ic42, &
560  ic51, ic52, ic61, ic62, ic71, ic72, ic81, ic82, &
561  dal1, dal2, dal3, af11, &
562  awg1, awg2, awg3, awg4, awg5, awg6, awg7, awg8, &
564  USE w3odatmd, ONLY: ndst, ndse
565 #ifdef W3_S
566  USE w3servmd, ONLY: strace
567 #endif
568  !/
569  IMPLICIT NONE
570  !/
571  !/ ------------------------------------------------------------------- /
572  !/ Parameter list
573  !/
574  INTEGER, INTENT(IN) :: IMOD
575  !/
576  !/ Local parameters
577  !/
578  INTEGER :: IFR, ITH, ISP, ITHP, ITHP1, ITHM, &
579  ITHM1,IFRP, IFRP1, IFRM, IFRM1
580  INTEGER, ALLOCATABLE :: IF1(:), IF2(:), IF3(:), IF4(:), &
581  IF5(:), IF6(:), IF7(:), IF8(:), &
582  IT1(:), IT2(:), IT3(:), IT4(:), &
583  IT5(:), IT6(:), IT7(:), IT8(:)
584 #ifdef W3_S
585  INTEGER, SAVE :: IENT = 0
586 #endif
587  REAL :: DELTH3, DELTH4, LAMM2, LAMP2, CTHP, &
588  WTHP, WTHP1, CTHM, WTHM, WTHM1, &
589  XFRLN, WFRP, WFRP1, WFRM, WFRM1, FR, &
590  AF11A
591  !/
592  !/ ------------------------------------------------------------------- /
593  !/
594 #ifdef W3_S
595  CALL strace (ient, 'INSNL1')
596 #endif
597 #ifdef W3_T
598  WRITE (ndst,9000) imod
599 #endif
600  !
601  nfr = nk
602  !
603  ! 1. Internal angles of quadruplet.
604  !
605  lamm2 = (1.-lam)**2
606  lamp2 = (1.+lam)**2
607  delth3 = acos( (lamm2**2+4.-lamp2**2) / (4.*lamm2) )
608  delth4 = asin(-sin(delth3)*lamm2/lamp2)
609  !
610  ! 2. Lambda dependend weight factors.
611  !
612  dal1 = 1. / (1.+lam)**4
613  dal2 = 1. / (1.-lam)**4
614  dal3 = 2. * dal1 * dal2
615  !
616  ! 3. Directional indices.
617  !
618  cthp = abs(delth4/dth)
619  ithp = int(cthp)
620  ithp1 = ithp + 1
621  wthp = cthp - real(ithp)
622  wthp1 = 1.- wthp
623  !
624  cthm = abs(delth3/dth)
625  ithm = int(cthm)
626  ithm1 = ithm + 1
627  wthm = cthm - real(ithm)
628  wthm1 = 1.- wthm
629  !
630  ! 4. Frequency indices.
631  !
632  xfrln = log(xfr)
633  !
634  ifrp = int( log(1.+lam) / xfrln )
635  ifrp1 = ifrp + 1
636  wfrp = (1.+lam - xfr**ifrp) / (xfr**ifrp1 - xfr**ifrp)
637  wfrp1 = 1. - wfrp
638  !
639  ifrm = int( log(1.-lam) / xfrln )
640  ifrm1 = ifrm - 1
641  wfrm = (xfr**ifrm -(1.-lam)) / (xfr**ifrm - xfr**ifrm1)
642  wfrm1 = 1. - wfrm
643  !
644  ! 5. Range of calculations
645  !
646  nfrhgh = nfr + ifrp1 - ifrm1
647  nfrchg = nfr - ifrm1
648  nspecy = nfrhgh * nth
649  nspecx = nfrchg * nth
650  !
651  ! 6. Allocate arrays or check array sizes
652  !
653  CALL w3dmnl ( imod, ndse, ndst, nspec, nspecx )
654  !
655  ALLOCATE ( if1(nfrchg), if2(nfrchg), if3(nfrchg), if4(nfrchg), &
656  if5(nfrchg), if6(nfrchg), if7(nfrchg), if8(nfrchg), &
657  it1(nth), it2(nth), it3(nth), it4(nth), &
658  it5(nth), it6(nth), it7(nth), it8(nth) )
659  !
660  ! 7. Spectral addresses
661  !
662  DO ifr=1, nfrchg
663  if1(ifr) = ifr+ifrp
664  if2(ifr) = ifr+ifrp1
665  if3(ifr) = max( 0 , ifr+ifrm )
666  if4(ifr) = max( 0 , ifr+ifrm1 )
667  if5(ifr) = max( 0 , ifr-ifrp )
668  if6(ifr) = max( 0 , ifr-ifrp1 )
669  if7(ifr) = ifr-ifrm
670  if8(ifr) = ifr-ifrm1
671  END DO
672  !
673  DO ith=1, nth
674  it1(ith) = ith + ithp
675  it2(ith) = ith + ithp1
676  it3(ith) = ith + ithm
677  it4(ith) = ith + ithm1
678  it5(ith) = ith - ithp
679  it6(ith) = ith - ithp1
680  it7(ith) = ith - ithm
681  it8(ith) = ith - ithm1
682  IF ( it1(ith).GT.nth) it1(ith) = it1(ith) - nth
683  IF ( it2(ith).GT.nth) it2(ith) = it2(ith) - nth
684  IF ( it3(ith).GT.nth) it3(ith) = it3(ith) - nth
685  IF ( it4(ith).GT.nth) it4(ith) = it4(ith) - nth
686  IF ( it5(ith).LT. 1 ) it5(ith) = it5(ith) + nth
687  IF ( it6(ith).LT. 1 ) it6(ith) = it6(ith) + nth
688  IF ( it7(ith).LT. 1 ) it7(ith) = it7(ith) + nth
689  IF ( it8(ith).LT. 1 ) it8(ith) = it8(ith) + nth
690  END DO
691  !
692  DO isp=1, nspecx
693  ifr = 1 + (isp-1)/nth
694  ith = 1 + mod(isp-1,nth)
695  ip11(isp) = it2(ith) + (if2(ifr)-1)*nth
696  ip12(isp) = it1(ith) + (if2(ifr)-1)*nth
697  ip13(isp) = it2(ith) + (if1(ifr)-1)*nth
698  ip14(isp) = it1(ith) + (if1(ifr)-1)*nth
699  im11(isp) = it8(ith) + (if4(ifr)-1)*nth
700  im12(isp) = it7(ith) + (if4(ifr)-1)*nth
701  im13(isp) = it8(ith) + (if3(ifr)-1)*nth
702  im14(isp) = it7(ith) + (if3(ifr)-1)*nth
703  ip21(isp) = it6(ith) + (if2(ifr)-1)*nth
704  ip22(isp) = it5(ith) + (if2(ifr)-1)*nth
705  ip23(isp) = it6(ith) + (if1(ifr)-1)*nth
706  ip24(isp) = it5(ith) + (if1(ifr)-1)*nth
707  im21(isp) = it4(ith) + (if4(ifr)-1)*nth
708  im22(isp) = it3(ith) + (if4(ifr)-1)*nth
709  im23(isp) = it4(ith) + (if3(ifr)-1)*nth
710  im24(isp) = it3(ith) + (if3(ifr)-1)*nth
711  END DO
712  !
713  DO isp=1, nspec
714  ifr = 1 + (isp-1)/nth
715  ith = 1 + mod(isp-1,nth)
716  ic11(isp) = it6(ith) + (if6(ifr)-1)*nth
717  ic21(isp) = it5(ith) + (if6(ifr)-1)*nth
718  ic31(isp) = it6(ith) + (if5(ifr)-1)*nth
719  ic41(isp) = it5(ith) + (if5(ifr)-1)*nth
720  ic51(isp) = it4(ith) + (if8(ifr)-1)*nth
721  ic61(isp) = it3(ith) + (if8(ifr)-1)*nth
722  ic71(isp) = it4(ith) + (if7(ifr)-1)*nth
723  ic81(isp) = it3(ith) + (if7(ifr)-1)*nth
724  ic12(isp) = it2(ith) + (if6(ifr)-1)*nth
725  ic22(isp) = it1(ith) + (if6(ifr)-1)*nth
726  ic32(isp) = it2(ith) + (if5(ifr)-1)*nth
727  ic42(isp) = it1(ith) + (if5(ifr)-1)*nth
728  ic52(isp) = it8(ith) + (if8(ifr)-1)*nth
729  ic62(isp) = it7(ith) + (if8(ifr)-1)*nth
730  ic72(isp) = it8(ith) + (if7(ifr)-1)*nth
731  ic82(isp) = it7(ith) + (if7(ifr)-1)*nth
732  END DO
733  !
734  DEALLOCATE ( if1, if2, if3, if4, if5, if6, if7, if8, &
735  it1, it2, it3, it4, it5, it6, it7, it8 )
736  !
737  ! 8. Fill scaling array (f**11)
738  !
739  DO ifr=1, nfr
740  af11a = (sig(ifr)*tpiinv)**11
741  DO ith=1, nth
742  af11(ith+(ifr-1)*nth) = af11a
743  END DO
744  END DO
745  !
746  fr = sig(nfr)*tpiinv
747  DO ifr=nfr+1, nfrchg
748  fr = fr * xfr
749  af11a = fr**11
750  DO ith=1, nth
751  af11(ith+(ifr-1)*nth) = af11a
752  END DO
753  END DO
754  !
755  ! 9. Interpolation weights
756  !
757  awg1 = wthp * wfrp
758  awg2 = wthp1 * wfrp
759  awg3 = wthp * wfrp1
760  awg4 = wthp1 * wfrp1
761  awg5 = wthm * wfrm
762  awg6 = wthm1 * wfrm
763  awg7 = wthm * wfrm1
764  awg8 = wthm1 * wfrm1
765  !
766  swg1 = awg1**2
767  swg2 = awg2**2
768  swg3 = awg3**2
769  swg4 = awg4**2
770  swg5 = awg5**2
771  swg6 = awg6**2
772  swg7 = awg7**2
773  swg8 = awg8**2
774  !
775  RETURN
776  !
777  ! Formats
778  !
779 #ifdef W3_T
780 9000 FORMAT (' TEST INSNL1 : IMOD :',i4)
781 #endif
782  !/
783  !/ End of INSNL1 ----------------------------------------------------- /
784  !/

References w3adatmd::af11, w3adatmd::awg1, w3adatmd::awg2, w3adatmd::awg3, w3adatmd::awg4, w3adatmd::awg5, w3adatmd::awg6, w3adatmd::awg7, w3adatmd::awg8, w3adatmd::dal1, w3adatmd::dal2, w3adatmd::dal3, w3gdatmd::dth, w3adatmd::ic11, w3adatmd::ic12, w3adatmd::ic21, w3adatmd::ic22, w3adatmd::ic31, w3adatmd::ic32, w3adatmd::ic41, w3adatmd::ic42, w3adatmd::ic51, w3adatmd::ic52, w3adatmd::ic61, w3adatmd::ic62, w3adatmd::ic71, w3adatmd::ic72, w3adatmd::ic81, w3adatmd::ic82, w3adatmd::im11, w3adatmd::im12, w3adatmd::im13, w3adatmd::im14, w3adatmd::im21, w3adatmd::im22, w3adatmd::im23, w3adatmd::im24, w3adatmd::ip11, w3adatmd::ip12, w3adatmd::ip13, w3adatmd::ip14, w3adatmd::ip21, w3adatmd::ip22, w3adatmd::ip23, w3adatmd::ip24, w3gdatmd::lam, w3odatmd::ndse, w3odatmd::ndst, w3adatmd::nfr, w3adatmd::nfrchg, w3adatmd::nfrhgh, w3gdatmd::nk, w3gdatmd::nspec, w3adatmd::nspecx, w3adatmd::nspecy, w3gdatmd::nth, w3gdatmd::sig, w3servmd::strace(), w3adatmd::swg1, w3adatmd::swg2, w3adatmd::swg3, w3adatmd::swg4, w3adatmd::swg5, w3adatmd::swg6, w3adatmd::swg7, w3adatmd::swg8, constants::tpiinv, w3adatmd::w3dmnl(), and w3gdatmd::xfr.

Referenced by w3iogrmd::w3iogr().

◆ insnlgqm()

subroutine w3snl1md::insnlgqm

Definition at line 1500 of file w3snl1md.F90.

1500  !/
1501  !/ +-----------------------------------+
1502  !/ | WAVEWATCH III NOAA/NCEP |
1503  !/ | E. Gagnaire-Renou & |
1504  !/ | M. Benoit |
1505  !/ | S. Mostafa Siadatamousavi |
1506  !/ | M. Beyramzadeh |
1507  !/ | FORTRAN 90 |
1508  !/ | Last update : 20-Nov-2022 |
1509  !/ +-----------------------------------+
1510  !/
1511  !/ 20-Nov-2022 : Merging with NL2 in WW3. ( version 7.00 )
1512  !/
1513  ! 1. Purpose :
1514  !
1515  ! Preprocessing for nonlinear interactions (Xnl).
1516  !
1517  ! 2. Method :
1518  !
1519  ! See Xnl documentation.
1520  !
1521  ! 3. Parameters :
1522  !
1523  ! 4. Subroutines used :
1524  !
1525  ! Name Type Module Description
1526  ! ----------------------------------------------------------------
1527  ! STRACE Subr. W3SERVMD Subroutine tracing.
1528  ! Subr. GAULEG Gauss-Legendre weights
1529  ! xnl_init Subr. m_constants Xnl initialization routine.
1530  ! ----------------------------------------------------------------
1531  !
1532  ! 5. Called by :
1533  !
1534  ! Name Type Module Description
1535  ! ----------------------------------------------------------------
1536  ! W3IOGR Subr. W3IOGRMD Model definition file management.
1537  ! ----------------------------------------------------------------
1538  !
1539  ! 6. Error messages :
1540  !
1541  ! 7. Remarks :
1542  !
1543  ! 8. Structure :
1544  !
1545  ! - See source code.
1546  !
1547  ! 9. Switches :
1548  !
1549  ! !/S Enable subroutine tracing.
1550  !
1551  ! 10. Source code :
1552  !
1553  !/ ------------------------------------------------------------------- /
1554  USE constants, ONLY: grav
1555  USE w3gdatmd, ONLY: nk , nth , xfr , fr1, gqnf1, gqnt1, gqnq_om2, nltail, gqthrcou
1556 
1557 #ifdef W3_S
1558  CALL strace (ient, 'INSNLGQM')
1559 #endif
1560  IMPLICIT NONE
1561  !.....LOCAL VARIABLES
1562  INTEGER JF , JT , JF1 , JT1 , NF1P1 , IAUX , NT , NF , IK
1563  INTEGER IQ_TE1 , IQ_OM2 , LBUF , DIMBUF , IQ_OM1 , NQ_TE1 , NCONFM
1564 
1565  DOUBLE PRECISION EPSI_A, AUX , CCC , DENO , AAA , DP2SG , TAILF
1566  DOUBLE PRECISION V1 , V1_4 , DV1 , DTETAR , ELIM , RAISF
1567  DOUBLE PRECISION V2 , V2_4 , V3 , V3_4
1568  DOUBLE PRECISION W2 , W2_M , W2_1 , W_MIL , W_RAD
1569  DOUBLE PRECISION RK0 , XK0 , YK0 , RK1 , XK1 , YK1
1570  DOUBLE PRECISION RK2 , XK2P , YK2P , XK2M , YK2M
1571  DOUBLE PRECISION RK3 , XK3P , YK3P , XK3M , YK3M
1572  DOUBLE PRECISION D01P , C_D01P, S_D01P, D0AP , C_D0AP, S_D0AP
1573  DOUBLE PRECISION GA2P , C_GA2P, S_GA2P, GA3P , C_GA3P, S_GA3P, TWOPI, PI, SEUIL1 , SEUIL2 , SEUIL
1574  !
1575  !.....Variables related to the Gaussian quadratures
1576  DOUBLE PRECISION W_CHE_TE1, W_CHE_OM2, C_LEG_OM2
1577  !
1578  !.....Variables related to the configuration selection
1579  DOUBLE PRECISION TEST1 , TEST2
1580  DOUBLE PRECISION :: FREQ(NK)
1581  DOUBLE PRECISION, ALLOCATABLE :: F1SF(:) , X_CHE_TE1(:) , X_CHE_OM2(:) , X_LEG_OM2(:) , W_LEG_OM2(:) &
1582  , MAXCLA(:)
1583 
1584  pi = acos(-1.)
1585  lbuf = 500
1586  dimbuf = 2*lbuf+200
1587  twopi = 2.*pi
1588  !
1589  ! Defines some threshold values for filtering (See Gagnaire-Renou Thesis, p 52)
1590  !
1591  seuil1 = 1e10
1592  seuil2 = gqthrcou
1593 
1594  IF(gqnf1.EQ.14) iq_om1=1
1595  IF(gqnf1.EQ.26) iq_om1=2
1596  IF(gqnf1.EQ.11) iq_om1=3
1597  IF(gqnf1.EQ.40) iq_om1=4
1598  IF(gqnf1.EQ.11) iq_om1=3
1599  IF(gqnf1.EQ.40) iq_om1=4
1600  IF(gqnf1.EQ.20) iq_om1=7
1601  !
1602  ! Note by FA: not sure what the 5 and 6 cases correspond to
1603  !
1604  nq_te1 = gqnt1/2
1605  nconfm = gqnf1*gqnt1*gqnq_om2
1606 
1607  raisf = xfr
1608  nt = nth
1609  nf = nk
1610  dtetar = twopi/dble(nt)
1611 
1612  DO ik = 1,nk
1613  freq(ik) = fr1*raisf**(ik-1)
1614  ENDDO
1615 
1616  tailf = -nltail
1617 
1618  !===============ALLOCATE MATRICES=============================================
1619  if (Allocated(k_if2) ) then
1620  deallocate(k_if2)
1621  endif
1622  ALLOCATE(k_if2(gqnq_om2,gqnt1,gqnf1))
1623 
1624  if (Allocated(k_if3) ) then
1625  deallocate(k_if3)
1626  endif
1627  ALLOCATE(k_if3(gqnq_om2,gqnt1,gqnf1))
1628 
1629  if (Allocated(k_1p2p) ) then
1630  deallocate(k_1p2p)
1631  endif
1632  ALLOCATE(k_1p2p(gqnq_om2,gqnt1,gqnf1))
1633 
1634  if (Allocated(k_1p3m) ) then
1635  deallocate(k_1p3m)
1636  endif
1637  ALLOCATE(k_1p3m(gqnq_om2,gqnt1,gqnf1))
1638 
1639  if (Allocated(k_1p2m) ) then
1640  deallocate(k_1p2m)
1641  endif
1642  ALLOCATE(k_1p2m(gqnq_om2,gqnt1,gqnf1))
1643 
1644  if (Allocated(k_1p3p) ) then
1645  deallocate(k_1p3p)
1646  endif
1647  ALLOCATE(k_1p3p(gqnq_om2,gqnt1,gqnf1))
1648 
1649  if (Allocated(k_1m2p) ) then
1650  deallocate(k_1m2p)
1651  endif
1652  ALLOCATE(k_1m2p(gqnq_om2,gqnt1,gqnf1))
1653 
1654  if (Allocated(k_1m3m) ) then
1655  deallocate(k_1m3m)
1656  endif
1657  ALLOCATE(k_1m3m(gqnq_om2,gqnt1,gqnf1))
1658 
1659  if (Allocated(k_1m2m) ) then
1660  deallocate(k_1m2m)
1661  endif
1662  ALLOCATE(k_1m2m(gqnq_om2,gqnt1,gqnf1))
1663 
1664  if (Allocated(k_1m3p) ) then
1665  deallocate(k_1m3p)
1666  endif
1667  ALLOCATE(k_1m3p(gqnq_om2,gqnt1,gqnf1))
1668 
1669  if (Allocated(tb_v24) ) then
1670  deallocate(tb_v24)
1671  endif
1672  ALLOCATE(tb_v24(gqnq_om2,gqnt1,gqnf1))
1673 
1674  if (Allocated(tb_v34) ) then
1675  deallocate(tb_v34)
1676  endif
1677  ALLOCATE(tb_v34(gqnq_om2,gqnt1,gqnf1))
1678 
1679  if (Allocated(tb_tpm) ) then
1680  deallocate(tb_tpm)
1681  endif
1682  ALLOCATE(tb_tpm(gqnq_om2,gqnt1,gqnf1))
1683 
1684  if (Allocated(tb_tmp) ) then
1685  deallocate(tb_tmp)
1686  endif
1687  ALLOCATE(tb_tmp(gqnq_om2,gqnt1,gqnf1))
1688 
1689  if (Allocated(tb_fac) ) then
1690  deallocate(tb_fac)
1691  endif
1692  ALLOCATE(tb_fac(gqnq_om2,gqnt1,gqnf1))
1693 
1694  if (Allocated(k_if1) ) then
1695  deallocate(k_if1)
1696  endif
1697  ALLOCATE(k_if1(gqnf1))
1698 
1699  if (Allocated(k_1p) ) then
1700  deallocate(k_1p)
1701  endif
1702  ALLOCATE(k_1p(gqnt1,gqnf1))
1703 
1704  if (Allocated(k_1m) ) then
1705  deallocate(k_1m)
1706  endif
1707  ALLOCATE(k_1m(gqnt1,gqnf1))
1708 
1709  if (Allocated(tb_v14) ) then
1710  deallocate(tb_v14)
1711  endif
1712  ALLOCATE(tb_v14(gqnf1))
1713 
1714  if (Allocated(idconf) ) then
1715  deallocate(idconf)
1716  endif
1717  ALLOCATE(idconf(nconfm,3))
1718 
1719  !=======================================================================
1720  ! INITIALISATION OF AUXILIAIRY TABLES FOR SPECTRUM INTERPOLATION
1721  !=======================================================================
1722  if (Allocated(f_poin) ) then
1723  deallocate(f_poin)
1724  endif
1725  ALLOCATE(f_poin(dimbuf))
1726 
1727  if (Allocated(t_poin) ) then
1728  deallocate(t_poin)
1729  endif
1730  ALLOCATE(t_poin(dimbuf))
1731 
1732  if (Allocated(f_coef) ) then
1733  deallocate(f_coef)
1734  endif
1735  ALLOCATE(f_coef(dimbuf))
1736 
1737  if (Allocated(f_proj) ) then
1738  deallocate(f_proj)
1739  endif
1740  ALLOCATE(f_proj(dimbuf))
1741 
1742  if (Allocated(tb_sca) ) then
1743  deallocate(tb_sca)
1744  endif
1745  ALLOCATE(tb_sca(dimbuf))
1746 
1747 
1748  f_poin(:)=0
1749  t_poin(:)=0
1750  f_coef(:)=0.d0
1751  f_proj(:)=0.d0
1752  tb_sca(:)=0.0d0
1753 
1754  DO jf=1,lbuf
1755  f_poin(jf)=1
1756  f_coef(jf)=0.0d0
1757  f_proj(jf)=0.0d0
1758  ENDDO
1759  DO jf=1,nf
1760  iaux=lbuf+jf
1761  f_poin(iaux)=jf
1762  f_coef(iaux)=1.0d0
1763  f_proj(iaux)=1.0d0
1764  ENDDO
1765  aux=1.d0/raisf**tailf
1766  DO jf=1,lbuf
1767  iaux=lbuf+nf+jf
1768  f_poin(iaux)=nf
1769  f_coef(iaux)=aux**jf
1770  f_proj(iaux)=0.0d0
1771  ENDDO
1772  !
1773  DO jt=lbuf,1,-1
1774  t_poin(jt)=nt-mod(lbuf-jt,nt)
1775  ENDDO
1776  DO jt=1,nt
1777  t_poin(lbuf+jt)=jt
1778  ENDDO
1779  DO jt=1,lbuf
1780  t_poin(lbuf+nt+jt)=mod(jt-1,nt)+1
1781  ENDDO
1782  !======================================================================
1783  !
1784  !=======================================================================
1785  ! COMPUTES SCALE COEFFICIENTS FOR THE COUPLING COEFFICIENT
1786  ! Would be easier to pass these on from W3SRCE ???
1787  !=======================================================================
1788  dp2sg=twopi*twopi/grav
1789  DO jf=1,lbuf
1790  aux=freq(1)/raisf**(lbuf-jf+1)
1791  tb_sca(jf)=(dp2sg*aux**2)**6/(twopi**3*aux)
1792  ENDDO
1793  DO jf=1,nf
1794  tb_sca(lbuf+jf)=(dp2sg*freq(jf)**2)**6/(twopi**3*freq(jf))
1795  ENDDO
1796  DO jf=1,lbuf
1797  iaux=lbuf+nf+jf
1798  aux=freq(nf)*raisf**jf
1799  tb_sca(iaux)=(dp2sg*aux**2)**6/(twopi**3*aux)
1800  ENDDO
1801  !=======================================================================
1802  !
1803  !=======================================================================
1804  ! COMPUTES VALUES FOR GAUSSIAN QUADRATURES
1805  !=======================================================================
1806  if (Allocated(x_che_te1) ) then
1807  deallocate(x_che_te1)
1808  endif
1809  ALLOCATE(x_che_te1(1:nq_te1),x_che_om2(1:gqnq_om2))
1810 
1811  if (Allocated(x_leg_om2) ) then
1812  deallocate(x_leg_om2)
1813  endif
1814  ALLOCATE(x_leg_om2(1:gqnq_om2),w_leg_om2(1:gqnq_om2))
1815  !
1816  !.....Abscissa and weight (constant) for Gauss-Chebyshev
1817  DO iq_te1=1,nq_te1
1818  x_che_te1(iq_te1)=cos(pi*(dble(iq_te1)-0.5d0)/dble(nq_te1))
1819  ENDDO
1820  w_che_te1=pi/dble(nq_te1)
1821  DO iq_om2=1,gqnq_om2
1822  x_che_om2(iq_om2)=cos(pi*(dble(iq_om2)-0.5d0)/dble(gqnq_om2))
1823  ENDDO
1824  w_che_om2=pi/dble(gqnq_om2)
1825  !
1826  !.....Abscissa et weight for Gauss-Legendre
1827  CALL gauleg( w_leg_om2 , x_leg_om2 , gqnq_om2 )
1828  DO iq_om2=1,gqnq_om2
1829  x_leg_om2(iq_om2)=0.25d0*(1.d0+x_leg_om2(iq_om2))**2
1830  ENDDO
1831  !=======================================================================
1832  !
1833  !
1834  !=======================================================================
1835  ! COMPUTES VALUES OF RATIO F1/F AS FUNCTION OF THE IQ_OM1 INDICATOR
1836  !=======================================================================
1837  nf1p1=gqnf1+1
1838  if (Allocated(f1sf) ) then
1839  deallocate(f1sf)
1840  endif
1841  ALLOCATE(f1sf(1:nf1p1))
1842 
1843  CALL f1f1f1 ( f1sf , gqnf1 , iq_om1)
1844  !=======================================================================
1845  !
1846  ! ==================================================
1847  ! STARTS LOOP 1 OVER THE RATIOS F1/F0
1848  ! ==================================================
1849  DO jf1=1,gqnf1
1850  ! ---------Computes and stores v1=f1/f0 and v1**4
1851  v1=(f1sf(jf1+1)+f1sf(jf1))/2.d0
1852  k_if1(jf1)=nint(dble(lbuf)+log(v1)/log(raisf))
1853  v1_4=v1**4
1854  tb_v14(jf1)=v1_4
1855  ! ---------Computes and stores dv1=df1/f0
1856  dv1=f1sf(jf1+1)-f1sf(jf1)
1857  ! ---------Computes the A parameter
1858  aaa=((1.d0+v1)**4-4.d0*(1.d0+v1_4))/(8.d0*v1**2)
1859  !
1860  ! =================================================
1861  ! STARTS LOOP 2 OVER THE DELTA_1+ VALUES
1862  ! =================================================
1863  DO jt1=1,gqnt1
1864  !
1865  !......Computes the Delta1+ values (=Theta_1-Theta_0) between 0 and Pi.
1866  IF (jt1.LE.nq_te1) THEN
1867  ! ---------First interval : X from -1 to A
1868  iq_te1=jt1
1869  c_d01p=(-1.d0+aaa)/2.d0+(1.d0+aaa)/2.d0*x_che_te1(iq_te1)
1870  ccc=dv1*sqrt((aaa-c_d01p)/(1.d0-c_d01p))*w_che_te1
1871  ELSE
1872  ! ---------Second interval : X from A to 1
1873  iq_te1=jt1-nq_te1
1874  c_d01p=( 1.d0+aaa)/2.d0+(1.d0-aaa)/2.d0*x_che_te1(iq_te1)
1875  ccc=dv1*sqrt((c_d01p-aaa)/(1.d0+c_d01p))*w_che_te1
1876  ENDIF
1877  s_d01p=sqrt(1.d0-c_d01p*c_d01p)
1878  d01p =acos(c_d01p)
1879  k_1p(jt1,jf1)=lbuf+nint(d01p/dtetar)
1880  k_1m(jt1,jf1)=lbuf-nint(d01p/dtetar)
1881  !
1882  ! ---------Computes Epsilon_a
1883  epsi_a=2.d0*sqrt(1.d0+v1_4+2.d0*v1*v1*c_d01p)/(1.d0+v1)**2
1884  ! ---------Computes Delta_A+ and its cosinus
1885  c_d0ap=(1.d0-v1_4+0.25d0*epsi_a**2*(1.d0+v1)**4) &
1886  /(epsi_a*(1.d0+v1)**2)
1887  s_d0ap=sqrt(1.0d0-c_d0ap*c_d0ap)
1888  d0ap = acos(c_d0ap)
1889  !
1890  !.......Integration over OMEGA2 depending on EPS_A
1891  IF (epsi_a.LT.1.d0) THEN
1892  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1893  !........Case of a single singularity (in OMEGA2-)
1894  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1895  w2_m=0.5d0*(1.d0-epsi_a/2.d0)
1896  w2_1=0.5d0
1897  !
1898  w_rad=w2_1-w2_m
1899  c_leg_om2=sqrt(w_rad)
1900  !
1901  ! ----------------------------------------------------
1902  !........STARTS LOOP 3 OVER OMEGA_2 (CASE Epsilon_A < 1)
1903  !........Case of a single singularity (in OMEGA2-)
1904  !........Integration over OMEGA2 via GAUSS-LEGENDRE quadrature
1905  ! ----------------------------------------------------
1906  DO iq_om2=1,gqnq_om2
1907  ! ---------Computes W2, V2, and V3
1908  w2=w2_m+w_rad*x_leg_om2(iq_om2)
1909  v2=w2*(1.d0+v1)
1910  v2_4=v2**4
1911  tb_v24(iq_om2,jt1,jf1)=v2_4
1912  k_if2(iq_om2,jt1,jf1) = nint(dble(lbuf) &
1913  + log(v2)/log(raisf))
1914  v3=1.d0+v1-v2
1915  v3_4=v3**4
1916  tb_v34(iq_om2,jt1,jf1)=v3_4
1917  k_if3(iq_om2,jt1,jf1) = nint(dble(lbuf) &
1918  + log(v3)/log(raisf))
1919  ! ---------Computes Gamma_2+ et Gamma_3+ angles
1920  c_ga2p=(epsi_a**2/4.d0+w2**4-(1.d0-w2)**4)/(epsi_a*w2*w2)
1921  c_ga2p=max(min(c_ga2p,1.d0),-1.d0)
1922  s_ga2p=sqrt(1.d0-c_ga2p*c_ga2p)
1923  ga2p =acos(c_ga2p)
1924  c_ga3p=(epsi_a**2/4.d0-w2**4+(1.d0-w2)**4)/epsi_a &
1925  /(1.d0-w2)**2
1926  c_ga3p=max(min(c_ga3p,1.d0),-1.d0)
1927  s_ga3p=sqrt(1.d0-c_ga3p*c_ga3p)
1928  ga3p =acos(c_ga3p)
1929  ! Shifting of the direction indexes - Config. +Delta1 (SIG=1)
1930  k_1p2p(iq_om2,jt1,jf1)=nint(( d0ap+ga2p)/dtetar &
1931  +dble(lbuf))
1932  k_1p3m(iq_om2,jt1,jf1)=nint(( d0ap-ga3p)/dtetar &
1933  +dble(lbuf))
1934  k_1p2m(iq_om2,jt1,jf1)=nint(( d0ap-ga2p)/dtetar &
1935  +dble(lbuf))
1936  k_1p3p(iq_om2,jt1,jf1)=nint(( d0ap+ga3p)/dtetar &
1937  +dble(lbuf))
1938  ! Shifting of the direction indexes - Config. -Delta1 (SIG=-1)
1939  k_1m2p(iq_om2,jt1,jf1)=nint((-d0ap+ga2p)/dtetar &
1940  +dble(lbuf))
1941  k_1m3m(iq_om2,jt1,jf1)=nint((-d0ap-ga3p)/dtetar &
1942  +dble(lbuf))
1943  k_1m2m(iq_om2,jt1,jf1)=nint((-d0ap-ga2p)/dtetar &
1944  +dble(lbuf))
1945  k_1m3p(iq_om2,jt1,jf1)=nint((-d0ap+ga3p)/dtetar &
1946  +dble(lbuf))
1947  !
1948  !.........Computes the coupling coefficients (only for Delta_1+ )
1949  rk0=1.d0
1950  rk1=v1*v1
1951  rk2=v2*v2
1952  rk3=(1.d0+v1-v2)**2
1953  xk0 = rk0
1954  yk0 = 0.0d0
1955  xk1 = rk1*c_d01p
1956  yk1 = rk1*s_d01p
1957  xk2p = rk2*(c_d0ap*c_ga2p-s_d0ap*s_ga2p)
1958  yk2p = rk2*(s_d0ap*c_ga2p+c_d0ap*s_ga2p)
1959  xk2m = rk2*(c_d0ap*c_ga2p+s_d0ap*s_ga2p)
1960  yk2m = rk2*(s_d0ap*c_ga2p-c_d0ap*s_ga2p)
1961  xk3p = rk3*(c_d0ap*c_ga3p-s_d0ap*s_ga3p)
1962  yk3p = rk3*(s_d0ap*c_ga3p+c_d0ap*s_ga3p)
1963  xk3m = rk3*(c_d0ap*c_ga3p+s_d0ap*s_ga3p)
1964  yk3m = rk3*(s_d0ap*c_ga3p-c_d0ap*s_ga3p)
1965  tb_tpm(iq_om2,jt1,jf1)=couple( xk0 , yk0 , xk1 , yk1 , xk2p , yk2p , xk3m , yk3m)
1966  tb_tmp(iq_om2,jt1,jf1)=couple( xk0 , yk0 , xk1 , yk1 , xk2m , yk2m , xk3p , yk3p)
1967  !
1968  !.........Computes the multiplicative coefficient for QNL4
1969  deno=2.d0*sqrt( (0.5d0*(1.d0+epsi_a/2.d0)-w2) &
1970  *((w2-0.5d0)**2+0.25d0*(1.d0+epsi_a)) &
1971  *((w2-0.5d0)**2+0.25d0*(1.d0-epsi_a)) )
1972  tb_fac(iq_om2,jt1,jf1)=1.d0/(deno*v1*w2*(1.d0-w2)) &
1973  /(1.d0+v1)**5 * w_leg_om2(iq_om2)*c_leg_om2* ccc
1974  ENDDO
1975  ! -----------------------------------------------
1976  !........END OF THE LOOP 3 OVER OMEGA_2 (CASE Epsilon_A < 1)
1977  ! -----------------------------------------------
1978  !
1979  ELSE
1980  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1981  !........STARTS LOOP 3 OVER OMEGA_2 (CASE Epsilon_A > 1)
1982  !........Case of two singularities (in OMEGA2- and OMEGA2_1)
1983  !........Integration over OMEGA2 via GAUSS-CHEBYSCHEV quadrature
1984  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1985  w2_m=0.5d0*(1.d0-epsi_a/2.d0)
1986  w2_1=0.5d0*(1.d0-sqrt(epsi_a-1.d0))
1987  !
1988  w_mil=(w2_m+w2_1)/2.d0
1989  w_rad=(w2_1-w2_m)/2.d0
1990  !
1991  DO iq_om2=1,gqnq_om2
1992  ! ---------Computes W2, V2, and V3
1993  w2=w_mil+w_rad*x_che_om2(iq_om2)
1994  v2=w2*(1.d0+v1)
1995  v2_4=v2**4
1996  tb_v24(iq_om2,jt1,jf1)=v2_4
1997  k_if2(iq_om2,jt1,jf1)=nint(dble(lbuf) &
1998  +log(v2)/log(raisf))
1999  v3=1.d0+v1-v2
2000  v3_4=v3**4
2001  tb_v34(iq_om2,jt1,jf1)=v3_4
2002  k_if3(iq_om2,jt1,jf1)=nint(dble(lbuf) &
2003  +log(v3)/log(raisf))
2004  ! ---------Computes Gamma_2+ et Gamma_3+ angles
2005  c_ga2p=(epsi_a**2/4.d0+w2**4-(1.d0-w2)**4)/(epsi_a*w2*w2)
2006  c_ga2p=max(min(c_ga2p,1.d0),-1.d0)
2007  s_ga2p=sqrt(1.d0-c_ga2p*c_ga2p)
2008  ga2p =acos(c_ga2p)
2009  c_ga3p=(epsi_a**2/4.d0-w2**4+(1.d0-w2)**4)/epsi_a &
2010  /(1.d0-w2)**2
2011  c_ga3p=max(min(c_ga3p,1.d0),-1.d0)
2012  s_ga3p=sqrt(1.d0-c_ga3p*c_ga3p)
2013  ga3p =acos(c_ga3p)
2014  ! Shifts the direction indexes - Config. +Delta1 (SIG=1)
2015  k_1p2p(iq_om2,jt1,jf1)=nint(( d0ap+ga2p)/dtetar &
2016  +dble(lbuf))
2017  k_1p3m(iq_om2,jt1,jf1)=nint(( d0ap-ga3p)/dtetar &
2018  +dble(lbuf))
2019  k_1p2m(iq_om2,jt1,jf1)=nint(( d0ap-ga2p)/dtetar &
2020  +dble(lbuf))
2021  k_1p3p(iq_om2,jt1,jf1)=nint(( d0ap+ga3p)/dtetar &
2022  +dble(lbuf))
2023  ! Shifts the direction indexes - Config. -Delta1 (SIG=-1)
2024  k_1m2p(iq_om2,jt1,jf1)=nint((-d0ap+ga2p)/dtetar &
2025  +dble(lbuf))
2026  k_1m3m(iq_om2,jt1,jf1)=nint((-d0ap-ga3p)/dtetar &
2027  +dble(lbuf))
2028  k_1m2m(iq_om2,jt1,jf1)=nint((-d0ap-ga2p)/dtetar &
2029  +dble(lbuf))
2030  k_1m3p(iq_om2,jt1,jf1)=nint((-d0ap+ga3p)/dtetar &
2031  +dble(lbuf))
2032  !
2033  !.........Computes the coupling coefficients (only for Delta_1+ )
2034  rk0=1.d0
2035  rk1=v1*v1
2036  rk2=v2*v2
2037  rk3=(1.d0+v1-v2)**2
2038  xk0 = rk0
2039  yk0 = 0.0d0
2040  xk1 = rk1*c_d01p
2041  yk1 = rk1*s_d01p
2042  xk2p = rk2*(c_d0ap*c_ga2p-s_d0ap*s_ga2p)
2043  yk2p = rk2*(s_d0ap*c_ga2p+c_d0ap*s_ga2p)
2044  xk2m = rk2*(c_d0ap*c_ga2p+s_d0ap*s_ga2p)
2045  yk2m = rk2*(s_d0ap*c_ga2p-c_d0ap*s_ga2p)
2046  xk3p = rk3*(c_d0ap*c_ga3p-s_d0ap*s_ga3p)
2047  yk3p = rk3*(s_d0ap*c_ga3p+c_d0ap*s_ga3p)
2048  xk3m = rk3*(c_d0ap*c_ga3p+s_d0ap*s_ga3p)
2049  yk3m = rk3*(s_d0ap*c_ga3p-c_d0ap*s_ga3p)
2050  tb_tpm(iq_om2,jt1,jf1)=couple( xk0 , yk0 , xk1 , yk1 , xk2p , yk2p , xk3m , yk3m)
2051  tb_tmp(iq_om2,jt1,jf1)=couple( xk0 , yk0 , xk1 , yk1 , xk2m , yk2m , xk3p , yk3p)
2052  !
2053  !.........Computes the multiplicative coefficient for QNL4
2054  deno=2.d0*sqrt( (0.5d0*(1.d0+epsi_a/2.d0)-w2) &
2055  *((w2-0.5d0)**2+0.25d0*(1.d0+epsi_a)) &
2056  *(0.5d0*(1.d0+sqrt(epsi_a-1.d0))-w2) )
2057  tb_fac(iq_om2,jt1,jf1)=1.d0/(deno*v1*w2*(1.d0-w2)) &
2058  /(1.d0+v1)**5 * w_che_om2* ccc
2059  !
2060  ENDDO
2061  ! -----------------------------------------------
2062  !........END OF LOOP 3 OVER OMEGA_2 (CASE Epsilon_A > 1)
2063  ! -----------------------------------------------
2064  !
2065  ENDIF
2066  ENDDO
2067  ! =================================================
2068  ! END OF LOOP 2 OVER THE DELTA_1+ VALUES
2069  ! =================================================
2070  !
2071  ENDDO
2072  ! ==================================================
2073  ! END OF LOOP 1 OVER THE F1/F0 RATIOS
2074  ! ==================================================
2075  DEALLOCATE(f1sf)
2076  DEALLOCATE(x_che_te1)
2077  DEALLOCATE(x_che_om2)
2078  DEALLOCATE(x_leg_om2)
2079  DEALLOCATE(w_leg_om2)
2080 
2081  ! ===========================================================
2082  ! POST-PROCESSING TO ELIMINATE PART OF THE CONFIGURATIONS
2083  ! ===========================================================
2084  !
2085  !.....It looks, for every value of the ratio V1, for the maximum value
2086  !.....of FACTOR*COUPLING : it is stored in the local table NAXCLA(.)
2087  ! """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""""""
2088  ALLOCATE(maxcla(1:gqnf1))
2089  DO jf1=1,gqnf1
2090  aux=0.0d0
2091  DO jt1=1,gqnt1
2092  DO iq_om2=1,gqnq_om2
2093  aux=max(aux,tb_fac(iq_om2,jt1,jf1)*tb_tpm(iq_om2,jt1,jf1),tb_fac(iq_om2,jt1,jf1)*tb_tmp(iq_om2,jt1,jf1))
2094  ENDDO
2095  ENDDO
2096  maxcla(jf1)=aux
2097  ENDDO
2098  !
2099  !.....It looks for the max V1 value
2100  ! """"""""""""""""""""""""""""""""""""""""""""""""""""""""""""
2101  aux=0.0d0
2102  DO jf1=1,gqnf1
2103  IF (maxcla(jf1).GT.aux) aux=maxcla(jf1)
2104  ENDDO
2105 
2106  test1=seuil1*aux
2107  !
2108  !.....Set to zero the coupling coefficients not used
2109  ! """""""""""""""""""""""""""""""""""""""""""""""""""""
2110  nconf=0
2111  DO jf1=1,gqnf1
2112  test2 =seuil2*maxcla(jf1)
2113  DO jt1=1,gqnt1
2114  DO iq_om2=1,gqnq_om2
2115  aaa=tb_fac(iq_om2,jt1,jf1)*tb_tpm(iq_om2,jt1,jf1)
2116  ccc=tb_fac(iq_om2,jt1,jf1)*tb_tmp(iq_om2,jt1,jf1)
2117  IF ((aaa.GT.test1.OR.aaa.GT.test2).OR. &
2118  (ccc.GT.test1.OR.ccc.GT.test2)) THEN
2119  nconf=nconf+1
2120  idconf(nconf,1)=jf1
2121  idconf(nconf,2)=jt1
2122  idconf(nconf,3)=iq_om2
2123  ENDIF
2124 #ifdef W3_TGQM
2125  WRITE(993,*) nconf,jf1,jt1,iq_om2,aaa,ccc,(aaa.GT.test1.OR.aaa.GT.test2), &
2126  (ccc.GT.test1.OR.ccc.GT.test2)
2127 #endif
2128  ENDDO
2129  ENDDO
2130  ENDDO
2131  DEALLOCATE(maxcla)
2132  !
2133  !..... counts the fraction of the eliminated configurations
2134  elim=(1.d0-dble(nconf)/dble(nconfm))*100.d0
2135 #ifdef W3_TGQM
2136  WRITE(994,*) 'NCONF, ELIM FRACTION:',nconf,elim
2137 #endif

References couple(), f1f1f1(), f_coef, f_poin, f_proj, w3gdatmd::fr1, gauleg(), w3gdatmd::gqnf1, w3gdatmd::gqnq_om2, w3gdatmd::gqnt1, w3gdatmd::gqthrcou, constants::grav, idconf, k_1m, k_1m2m, k_1m2p, k_1m3m, k_1m3p, k_1p, k_1p2m, k_1p2p, k_1p3m, k_1p3p, k_if1, k_if2, k_if3, nconf, w3gdatmd::nk, w3gdatmd::nltail, w3gdatmd::nth, t_poin, tb_fac, tb_sca, tb_tmp, tb_tpm, tb_v14, tb_v24, tb_v34, and w3gdatmd::xfr.

Referenced by w3iogrmd::w3iogr().

◆ w3snl1()

subroutine w3snl1md::w3snl1 ( real, dimension(nspec), intent(in)  A,
real, dimension(nk), intent(in)  CG,
real, intent(in)  KDMEAN,
real, dimension(nspec), intent(out)  S,
real, dimension(nspec), intent(out)  D 
)

Calculate nonlinear interactions and the diagonal term of its derivative.

Parameters
[in]AAction spectrum A(ISP) as a function of direction (rad) and wavenumber.
[in]CGGroup velocities (dimension NK).
[in]KDMEANMean relative depth.
[out]SSource term.
[out]DDiagonal term of derivative.
Author
H. L. Tolman
Date
06-Jun-2018

Definition at line 115 of file w3snl1md.F90.

115  !/
116  !/ +-----------------------------------+
117  !/ | WAVEWATCH III NOAA/NCEP |
118  !/ | H. L. Tolman |
119  !/ | FORTRAN 90 |
120  !/ | Last update : 06-Jun-2018 |
121  !/ +-----------------------------------+
122  !/
123  !/ 12-Jun-1996 : Final FORTRAN 77 ( version 1.18 )
124  !/ 04-Feb-2000 : Upgrade to FORTRAN 90 ( version 2.00 )
125  !/ 09-May-2002 : Switch clean up. ( version 2.21 )
126  !/ 24-Dec-2004 : Multiple grid version. ( version 3.06 )
127  !/ 03-Sep-2012 : Clean up of test output T0, T1 ( version 4.07 )
128  !/ 06-Jun-2018 : Add optional DEBUGSRC ( version 6.04 )
129  !/
130  ! 1. Purpose :
131  !
132  ! Calculate nonlinear interactions and the diagonal term of
133  ! its derivative.
134  !
135  ! 2. Method :
136  !
137  ! Discrete interaction approximation. (Hasselmann and Hasselmann
138  ! 1985; WAMDI group 1988)
139  !
140  ! The DIA is applied to the energy spectrum (instead of the action
141  ! spectrum), for which is was originally developped. Because the
142  ! frequency grid is invariant, the nonlinear interactions are
143  ! calculated for the frequency spectrum, as in WAM. This requires
144  ! only a single set of interpolation data which can be applied
145  ! throughout the spatial domain. For deep water this is idenitical
146  ! to a direct application to the wavenumber spectrum, for shallow
147  ! water it is not. As the shallow water correction is nothing but
148  ! a crude approximation, the choice between spectra is expected to
149  ! be irrelevant.
150  !
151  ! The nonlinear interactions are calculated for two "mirror image"
152  ! quadruplets as described in the manual. The central bin of these
153  ! quadruples is placed on the discrete complonents of the spectrum,
154  ! which requires interpolation to obtain other eneregy densities.
155  ! The figure below defines the diferent basic counters and weights
156  ! necessary for this interpolation.
157  !
158  !
159  ! IFRM1 IFRM
160  ! 5 7 T |
161  ! ITHM1 +------+ H +
162  ! | | E | IFRP IFRP1
163  ! | \ | T | 3 1
164  ! ITHM +------+ A + +---------+ ITHP1
165  ! 6 \8 | | |
166  ! | | / |
167  ! \ + +---------+ ITHP
168  ! | /4 2
169  ! \ | /
170  ! -+-----+------+-------#--------+---------+----------+
171  ! / | \ FREQ.
172  ! | \4 2
173  ! / + +---------+ ITHP
174  ! | | \ |
175  ! 6 /8 | | |
176  ! ITHM +------+ + +---------+ ITHP1
177  ! | \ | | 3 1
178  ! | | | IFRP IFRP1
179  ! ITHM1 +------+ +
180  ! 5 7 |
181  !
182  ! To create long vector loops and to efficiently deal with the
183  ! closed nature of the directional space, the relative counters
184  ! above are replaced by complete addresses stored in 32 arrays
185  ! (see section 3 and INSNL1). The interaction are furthermore
186  ! calucated for an extended spectrum, making it unnecessary to
187  ! introduce extra weight factors for low and high frequencies.
188  ! Therefore low and high frequencies are added to the local
189  ! (auxiliary) spectrum as illustraed below.
190  !
191  !
192  ! ^ +---+---------------------+---------+- NTH
193  ! | | : : |
194  ! | : : |
195  ! d | 2 : original spectrum : 1 |
196  ! i | : : |
197  ! r | : : |
198  ! +---+---------------------+---------+- 1
199  ! Frequencies --> ^
200  ! IFR = 0 1 NFR | NFRHGH
201  ! |
202  ! NFRCHG
203  !
204  ! where : 1 : Extra tail added beyond NFR
205  ! 2 : Empty bins at low frequencies
206  !
207  ! NFRHGH = NFR + IFRP1 - IFRM1
208  ! NFRCHG = NFR - IFRM1
209  !
210  ! All counters and arrays are set in INSNL1. See also section 3
211  ! and section 8.
212  !
213  ! 3. Parameters :
214  !
215  ! Parameter list
216  ! ----------------------------------------------------------------
217  ! A R.A. I Action spectrum A(ISP) as a function of
218  ! direction (rad) and wavenumber.
219  ! CG R.A. I Group velocities (dimension NK).
220  ! KDMEAN Real I Mean relative depth.
221  ! S R.A. O Source term. *)
222  ! D R.A. O Diagonal term of derivative. *)
223  ! ----------------------------------------------------------------
224  ! *) 1-D array with dimension NTH*NK
225  !
226  ! 4. Subroutines used :
227  !
228  ! Name Type Module Description
229  ! ----------------------------------------------------------------
230  ! STRACE Subr. W3SERVMD Subroutine tracing.
231  ! PRT2DS Subr. W3ARRYMD Print plot of spectra.
232  ! OUTMAT Subr. W3WRRYMD Print out 2D matrix.
233  ! ----------------------------------------------------------------
234  !
235  ! 5. Called by :
236  !
237  ! Name Type Module Description
238  ! ----------------------------------------------------------------
239  ! W3SRCE Subr. W3SRCEMD Source term integration.
240  ! W3EXPO Subr. N/A Point output post-processor.
241  ! GXEXPO Subr. N/A GrADS point output post-processor.
242  ! ----------------------------------------------------------------
243  !
244  ! 6. Error messages :
245  !
246  ! None.
247  !
248  ! 7. Remarks :
249  !
250  ! None.
251  !
252  ! 8. Structure :
253  !
254  ! -------------------------------------------
255  ! 1. Calculate proportionality constant.
256  ! 2. Prepare auxiliary spectrum
257  ! 3. Calculate (unfolded) interactions
258  ! a Energy at interacting bins
259  ! b Contribution to interactions
260  ! c Fold interactions to side angles
261  ! 4. Put source and diagonal term together
262  ! -------------------------------------------
263  !
264  ! 9. Switches :
265  !
266  ! !/S Enable subroutine tracing.
267  ! !/T Enable general test output.
268  ! !/T0 2-D print plot of source term.
269  ! !/T1 Print arrays.
270  !
271  ! 10. Source code :
272  !
273  !/ ------------------------------------------------------------------- /
274  !/
275  USE constants
276  USE w3gdatmd, ONLY: nk, nth, nspec, sig, fachfe, &
278  USE w3adatmd, ONLY: nfr, nfrhgh, nfrchg, nspecx, nspecy, &
279  ip11, ip12, ip13, ip14, im11, im12, im13, im14, &
280  ip21, ip22, ip23, ip24, im21, im22, im23, im24, &
281  ic11, ic12, ic21, ic22, ic31, ic32, ic41, ic42, &
282  ic51, ic52, ic61, ic62, ic71, ic72, ic81, ic82, &
283  dal1, dal2, dal3, af11, &
284  awg1, awg2, awg3, awg4, awg5, awg6, awg7, awg8, &
286 #ifdef W3_T
287  USE w3odatmd, ONLY: ndst
288 #endif
289 #ifdef W3_T1
290  USE w3odatmd, ONLY: ndst
291 #endif
292 #ifdef W3_S
293  USE w3servmd, ONLY: strace
294 #endif
295 #ifdef W3_T0
296  USE w3arrymd, ONLY: prt2ds
297 #endif
298 #ifdef W3_T1
299  USE w3arrymd, ONLY: outmat
300 #endif
301  !
302  IMPLICIT NONE
303  !/
304  !/ ------------------------------------------------------------------- /
305  !/ Parameter list
306  !/
307  REAL, INTENT(IN) :: A(NSPEC), CG(NK), KDMEAN
308  REAL, INTENT(OUT) :: S(NSPEC), D(NSPEC)
309  !/
310  !/ ------------------------------------------------------------------- /
311  !/ Local parameters
312  !/
313  INTEGER :: ITH, IFR, ISP
314 #ifdef W3_S
315  INTEGER, SAVE :: IENT = 0
316 #endif
317  REAL :: X, X2, CONS, CONX, FACTOR, &
318  E00, EP1, EM1, EP2, EM2, &
319  SA1A, SA1B, SA2A, SA2B
320 #ifdef W3_T0
321  REAL :: SOUT(NK,NFR), DOUT(NK,NFR)
322 #endif
323  REAL :: UE (1-NTH:NSPECY), SA1 (1-NTH:NSPECX), &
324  SA2 (1-NTH:NSPECX), DA1C(1-NTH:NSPECX), &
325  DA1P(1-NTH:NSPECX), DA1M(1-NTH:NSPECX), &
326  DA2C(1-NTH:NSPECX), DA2P(1-NTH:NSPECX), &
327  DA2M(1-NTH:NSPECX), CON ( NSPEC )
328  !/
329  !/ ------------------------------------------------------------------- /
330  !/
331  ! initialisations
332  !
333 #ifdef W3_S
334  CALL strace (ient, 'W3SNL1')
335 #endif
336  !
337  ! 1. Calculate prop. constant --------------------------------------- *
338  !
339  x = max( kdcon*kdmean , kdmn )
340  x2 = max( -1.e15, snls3*x)
341  cons = snlc1 * ( 1. + snls1/x * (1.-snls2*x) * exp(x2) )
342  !
343 #ifdef W3_T
344  WRITE (ndst,9000) kdmean, cons
345 #endif
346  !
347  ! 2. Prepare auxiliary spectrum and arrays -------------------------- *
348  !
349  DO ifr=1, nfr
350  conx = tpiinv / sig(ifr) * cg(ifr)
351  DO ith=1, nth
352  isp = ith + (ifr-1)*nth
353  ue(isp) = a(isp) / conx
354  con(isp) = conx
355  END DO
356  END DO
357  !
358  DO ifr=nfr+1, nfrhgh
359  DO ith=1, nth
360  isp = ith + (ifr-1)*nth
361  ue(isp) = ue(isp-nth) * fachfe
362  END DO
363  END DO
364  !
365  DO isp=1-nth, 0
366  ue(isp) = 0.
367  sa1(isp) = 0.
368  sa2(isp) = 0.
369  da1c(isp) = 0.
370  da1p(isp) = 0.
371  da1m(isp) = 0.
372  da2c(isp) = 0.
373  da2p(isp) = 0.
374  da2m(isp) = 0.
375  END DO
376  !
377  ! 3. Calculate interactions for extended spectrum ------------------- *
378  !
379  DO isp=1, nspecx
380  !
381  ! 3.a Energy at interacting bins
382  !
383  e00 = ue(isp)
384  ep1 = awg1 * ue(ip11(isp)) + awg2 * ue(ip12(isp)) &
385  + awg3 * ue(ip13(isp)) + awg4 * ue(ip14(isp))
386  em1 = awg5 * ue(im11(isp)) + awg6 * ue(im12(isp)) &
387  + awg7 * ue(im13(isp)) + awg8 * ue(im14(isp))
388  ep2 = awg1 * ue(ip21(isp)) + awg2 * ue(ip22(isp)) &
389  + awg3 * ue(ip23(isp)) + awg4 * ue(ip24(isp))
390  em2 = awg5 * ue(im21(isp)) + awg6 * ue(im22(isp)) &
391  + awg7 * ue(im23(isp)) + awg8 * ue(im24(isp))
392  !
393  ! 3.b Contribution to interactions
394  !
395  factor = cons * af11(isp) * e00
396  !
397  sa1a = e00 * ( ep1*dal1 + em1*dal2 )
398  sa1b = sa1a - ep1*em1*dal3
399  sa2a = e00 * ( ep2*dal1 + em2*dal2 )
400  sa2b = sa2a - ep2*em2*dal3
401  !
402  sa1(isp) = factor * sa1b
403  sa2(isp) = factor * sa2b
404  !
405  da1c(isp) = cons * af11(isp) * ( sa1a + sa1b )
406  da1p(isp) = factor * ( dal1*e00 - dal3*em1 )
407  da1m(isp) = factor * ( dal2*e00 - dal3*ep1 )
408  !
409  da2c(isp) = cons * af11(isp) * ( sa2a + sa2b )
410  da2p(isp) = factor * ( dal1*e00 - dal3*em2 )
411  da2m(isp) = factor * ( dal2*e00 - dal3*ep2 )
412  !
413  END DO
414  !
415  ! 4. Put source and diagonal term together -------------------------- *
416  !
417  DO isp=1, nspec
418  !
419  s(isp) = con(isp) * ( - 2. * ( sa1(isp) + sa2(isp) ) &
420  + awg1 * ( sa1(ic11(isp)) + sa2(ic12(isp)) ) &
421  + awg2 * ( sa1(ic21(isp)) + sa2(ic22(isp)) ) &
422  + awg3 * ( sa1(ic31(isp)) + sa2(ic32(isp)) ) &
423  + awg4 * ( sa1(ic41(isp)) + sa2(ic42(isp)) ) &
424  + awg5 * ( sa1(ic51(isp)) + sa2(ic52(isp)) ) &
425  + awg6 * ( sa1(ic61(isp)) + sa2(ic62(isp)) ) &
426  + awg7 * ( sa1(ic71(isp)) + sa2(ic72(isp)) ) &
427  + awg8 * ( sa1(ic81(isp)) + sa2(ic82(isp)) ) )
428  !
429  d(isp) = - 2. * ( da1c(isp) + da2c(isp) ) &
430  + swg1 * ( da1p(ic11(isp)) + da2p(ic12(isp)) ) &
431  + swg2 * ( da1p(ic21(isp)) + da2p(ic22(isp)) ) &
432  + swg3 * ( da1p(ic31(isp)) + da2p(ic32(isp)) ) &
433  + swg4 * ( da1p(ic41(isp)) + da2p(ic42(isp)) ) &
434  + swg5 * ( da1m(ic51(isp)) + da2m(ic52(isp)) ) &
435  + swg6 * ( da1m(ic61(isp)) + da2m(ic62(isp)) ) &
436  + swg7 * ( da1m(ic71(isp)) + da2m(ic72(isp)) ) &
437  + swg8 * ( da1m(ic81(isp)) + da2m(ic82(isp)) )
438  !
439  END DO
440  !
441  ! ... Test output :
442  !
443 #ifdef W3_T0
444  DO ifr=1, nfr
445  DO ith=1, nth
446  isp = ith + (ifr-1)*nth
447  sout(ifr,ith) = s(isp) * tpi * sig(ifr) / cg(ifr)
448  dout(ifr,ith) = d(isp)
449  END DO
450  END DO
451  CALL prt2ds (ndst, nk, nk, nth, sout, sig(1:), ' ', 1., &
452  0.0, 0.001, 'Snl(f,t)', ' ', 'NONAME')
453  CALL prt2ds (ndst, nk, nk, nth, dout, sig(1:), ' ', 1., &
454  0.0, 0.001, 'Diag Snl', ' ', 'NONAME')
455 #endif
456  !
457 #ifdef W3_T1
458  CALL outmat (ndst, s, nth, nth, nk, 'Snl')
459  CALL outmat (ndst, d, nth, nth, nk, 'Diag Snl')
460 #endif
461  !
462  RETURN
463  !
464  ! Formats
465  !
466 #ifdef W3_T
467 9000 FORMAT (' TEST W3SNL1 : KDMEAN, CONS :',f8.2,f8.1)
468 #endif
469  !/
470  !/ End of W3SNL1 ----------------------------------------------------- /
471  !/

References w3adatmd::af11, w3adatmd::awg1, w3adatmd::awg2, w3adatmd::awg3, w3adatmd::awg4, w3adatmd::awg5, w3adatmd::awg6, w3adatmd::awg7, w3adatmd::awg8, w3adatmd::dal1, w3adatmd::dal2, w3adatmd::dal3, w3gdatmd::fachfe, w3adatmd::ic11, w3adatmd::ic12, w3adatmd::ic21, w3adatmd::ic22, w3adatmd::ic31, w3adatmd::ic32, w3adatmd::ic41, w3adatmd::ic42, w3adatmd::ic51, w3adatmd::ic52, w3adatmd::ic61, w3adatmd::ic62, w3adatmd::ic71, w3adatmd::ic72, w3adatmd::ic81, w3adatmd::ic82, w3adatmd::im11, w3adatmd::im12, w3adatmd::im13, w3adatmd::im14, w3adatmd::im21, w3adatmd::im22, w3adatmd::im23, w3adatmd::im24, w3adatmd::ip11, w3adatmd::ip12, w3adatmd::ip13, w3adatmd::ip14, w3adatmd::ip21, w3adatmd::ip22, w3adatmd::ip23, w3adatmd::ip24, w3gdatmd::kdcon, w3gdatmd::kdmn, w3odatmd::ndst, w3adatmd::nfr, w3adatmd::nfrchg, w3adatmd::nfrhgh, w3gdatmd::nk, w3gdatmd::nspec, w3adatmd::nspecx, w3adatmd::nspecy, w3gdatmd::nth, w3arrymd::outmat(), w3arrymd::prt2ds(), w3gdatmd::sig, w3gdatmd::snlc1, w3gdatmd::snls1, w3gdatmd::snls2, w3gdatmd::snls3, w3servmd::strace(), w3adatmd::swg1, w3adatmd::swg2, w3adatmd::swg3, w3adatmd::swg4, w3adatmd::swg5, w3adatmd::swg6, w3adatmd::swg7, w3adatmd::swg8, constants::tpi, and constants::tpiinv.

Referenced by gxexpo(), w3exnc(), w3expo(), and w3srcemd::w3srce().

◆ w3snlgqm()

subroutine w3snl1md::w3snlgqm ( real, dimension(nth,nk), intent(in)  A,
real, dimension(nk), intent(in)  CG,
real, dimension(nk), intent(in)  WN,
real, intent(in)  DEPTH,
real, dimension(nth,nk), intent(out)  TSTOTn,
real, dimension(nth,nk), intent(out)  TSDERn 
)

Definition at line 789 of file w3snl1md.F90.

789  ! This and the following routines are adapted to WW3 from TOMAWAC qnlin3.f
790  !***********************************************************************
791  ! TOMAWAC V6P1 24/06/2011
792  !***********************************************************************
793  !
794  !brief COMPUTES THE CONTRIBUTION OF THE NON-LINEAR INTERACTIONS
795  !+ SOURCE TERM BETWEEN QUADRUPLETS USING THE GQM METHOD
796  !+ ("GAUSSIAN QUADRATURE METHOD") PROPOSED BY LAVRENOV
797  !+ (2001)
798  !+
799  !+ PROCEDURE SPECIFIC TO THE CASE WHERE THE FREQUENCIES
800  !+ FOLLOW A GEOMETRICAL PROGRESSION AND THE DIRECTIONS
801  !+ ARE EVENLY DISTRIBUTED OVER [0;2.PI].
802  !
803  !note THIS SUBROUTINE USES THE OUTPUT FROM 'PRENL3' TO OPTIMISE
804  !+ THE COMPUTATIONS FOR DIA.
805  !
806  !reference LAVRENOV, I.V. (2001):
807  !+ "EFFECT OF WIND WAVE PARAMETER FLUCTUATION ON THE NONLINEAR
808  !+ SPECTRUM EVOLUTION". J. PHYS. OCEANOGR. 31, 861-873.
809  !
810  !history E. GAGNAIRE-RENOU
811  !+ 04/2011
812  !+ V6P1
813  !+ CREATED
814  !
815  !history G.MATTAROLO (EDF - LNHE)
816  !+ 24/06/2011
817  !+ V6P1
818  !+ Translation of French names of the variables in argument
819 
820  !
821  !/ Warning, contrary to the DIA routine, there is no extension to frequencies below IK=1
822  !/ as a result the first two frequencies are not fully treated.
823  !==================================================================================
824  ! This subroutine is same as qnlin3 in TOMWAC
825  USE constants, ONLY: tpi
826  USE w3gdatmd, ONLY: sig, nk , nth , dth, xfr, fr1, gqthrsat, gqamp
827 
828  IMPLICIT NONE
829 
830  REAL, intent(in) :: A(NTH,NK), CG(NK), WN(NK)
831  REAL, intent(in) :: DEPTH
832  REAL, intent(out) :: TSTOTn(NTH,NK), TSDERn(NTH,NK)
833 
834  INTEGER :: ITH,IK,NT,NF
835  REAL :: q_dfac, SATVAL(NK), SUME, ACCVAL, ACCMAX, AMPFAC
836  DOUBLE PRECISION :: RAISF, FREQ(NK)
837  DOUBLE PRECISION :: TSTOT(NTH,NK) , TSDER(NTH,NK), F(NTH,NK)
838  DOUBLE PRECISION :: TEMP
839 
840  !.....LOCAL VARIABLES
841  INTEGER JF , JT , JF1 , JT1 , IQ_OM2 &
842  , JFM0 , JFM1 , JFM2 , JFM3 , IXF1 , IXF2 &
843  , IXF3 , JFMIN , JFMAX , ICONF , LBUF
844  INTEGER KT1P , KT1M , JT1P , JT1M , KT1P2P, KT1P2M &
845  , KT1P3P, KT1P3M, KT1M2P, KT1M2M, KT1M3P, KT1M3M &
846  , JT1P2P, JT1P2M, JT1P3P, JT1P3M, JT1M2P, JT1M2M &
847  , JT1M3P, JT1M3M
848  DOUBLE PRECISION V1_4 , V2_4 , V3_4 , Q_2P3M, Q_2M3P, FACTOR &
849  , T_2P3M, T_2M3P, S_2P3M, S_2M3P, SCAL_T, T2P3M &
850  , T2M3P , SP0 , SP1P , SP1M , SP1P2P, SP1P2M &
851  , SP1P3P, SP1P3M, SP1M2P, SP1M2M, SP1M3P, SP1M3M &
852  , CF0 , CP0 , CF1 , CP1 , CF2 , CP2 &
853  , CF3 , CP3 , Q2PD0 , Q2PD1 , Q2PD2P, Q2PD3M &
854  , Q2MD0 , Q2MD1 , Q2MD2M, Q2MD3P ,AUX00 , AUX01 &
855  , AUX02 , AUX03 , AUX04 , AUX05 , SEUIL &
856  , AUX06 , AUX07 , AUX08 , AUX09 , AUX10 , FSEUIL
857 
858  nt = nth
859  nf = nk
860  lbuf = 500
861  seuil = 0.
862  raisf = xfr
863 
864  DO ik = 1,nk
865  freq(ik) = fr1*raisf**(ik-1)
866  ENDDO
867 
868  DO ith = 1,nth
869  DO ik = 1,nk
870  ! F is the E(f,theta) spectrum ...
871  f(ith,ik) = dble(a(ith,ik)*sig(ik))*dble(tpi)/dble(cg(ik))
872  ENDDO
873  ENDDO
874  ! CALL INSNLGQM
875  ! it returns: F_POIN , T_POIN , F_COEF , F_PROJ, TB_SCA , K_IF1, K_1P, k_1M , K_IF2
876  ! K_IF3, K_1P2P , K_1P3M , K_1P2M , K_1P3P , K_1M2P , K_1M3M , K_1M2M
877  ! K_1M3P , TB_V14 , TB_FAC , TB_V24 , TB_V34 , TB_TMP , TB_TPM , IDCONF, NCONF
878  !=======================================================================
879  ! COMPUTES THE GENERALIZED MIN AND MAX FREQUENCIES : INSTEAD OF GOING
880  ! FROM 1 TO NF IN FREQ(JF) FOR THE MAIN FREQUENCY, IT GOES FROM JFMIN
881  ! TO JFMAX
882  ! JFMIN IS GIVEN BY Fmin=FREQ(1) /Gamma_min
883  ! JFMAX IS GIVEN BY Fmax=FREQ(NF)*Gamma_max
884  ! TESTS HAVE SHOWN THAT IT CAN BE ASSUMED Gamma_min=1. (JFMIN=1) AND
885  ! Gamma_max=1.3 (JFMAX>NF) TO OBTAIN IMPROVED RESULTS
886  ! Note by Fabrice Ardhuin: this appears to give the difference in tail benaviour with Gerbrant's WRT
887  !=======================================================================
888  jfmin=max(1-int(log(1.0d0)/log(raisf)),1)
889  jfmax=min(nf+int(log(1.3d0)/log(raisf)),nk)
890  !
891  !=======================================================================
892  ! COMPUTES THE SPECTRUM THRESHOLD VALUES (BELOW WHICH QNL4 IS NOT
893  ! CALCULATED). THE THRESHOLD IS SET WITHIN 0 AND 1.
894  ! This was commented by FA
895  !=======================================================================
896  ! AUX00=0.0D0
897  ! DO JF=1,NF
898  ! DO JT=1,NT
899  ! IF (F(JT,JF).GT.AUX00) AUX00=F(JT,JF)
900  ! ENDDO
901  ! ENDDO
902  ! FSEUIL=AUX00*SEUIL
903 
904  tstot = 0.
905  tsder = 0.
906  !=======================================================================
907  accmax=0.
908  DO jf=jfmin,jfmax
909  sume=sum(f(:,jf))*dth
910  satval(jf) = sume*freq(jf)**5
911  accval = sume*freq(jf)**4
912  IF (accval.GT.accmax) accmax=accval
913  END DO
914 
915 
916  ! ==================================================
917  ! STARTS LOOP 1 OVER THE SELECTED CONFIGURATIONS
918  ! ==================================================
919  DO iconf=1,nconf
920  ! ---------selected configuration characteristics
921  jf1 =idconf(iconf,1)
922  jt1 =idconf(iconf,2)
923  iq_om2=idconf(iconf,3)
924  !
925  ! ---------Recovers V1**4=(f1/f0)**4
926  v1_4 =tb_v14(jf1)
927  ! ---------Recovers the shift of the frequency index on f1
928  ixf1 =k_if1(jf1)
929  ! ---------Recovers the direction indexes for Delat1
930  kt1p =k_1p(jt1,jf1)
931  kt1m =k_1m(jt1,jf1)
932  ! ---------Recovers V2**4=(f2/f0)**4 and V3**4=(f3/f0)**4
933  v2_4 =tb_v24(iq_om2,jt1,jf1)
934  v3_4 =tb_v34(iq_om2,jt1,jf1)
935  ! ---------Recovers the frequency indexes shift on f2 and f3
936  ixf2 =k_if2(iq_om2,jt1,jf1)
937  ixf3 =k_if3(iq_om2,jt1,jf1)
938  ! ---------Recovers the direction indexes shift
939  kt1p2p=k_1p2p(iq_om2,jt1,jf1)
940  kt1p2m=k_1p2m(iq_om2,jt1,jf1)
941  kt1p3p=k_1p3p(iq_om2,jt1,jf1)
942  kt1p3m=k_1p3m(iq_om2,jt1,jf1)
943  kt1m2p=k_1m2p(iq_om2,jt1,jf1)
944  kt1m2m=k_1m2m(iq_om2,jt1,jf1)
945  kt1m3p=k_1m3p(iq_om2,jt1,jf1)
946  kt1m3m=k_1m3m(iq_om2,jt1,jf1)
947  ! ---------Recovers the coupling coefficients
948  t2p3m =tb_tpm(iq_om2,jt1,jf1)
949  t2m3p =tb_tmp(iq_om2,jt1,jf1)
950  ! ---------Recovers the multiplicative factor of QNL4
951  factor=tb_fac(iq_om2,jt1,jf1)
952 
953  ! = = = = = = = = = = = = = = = = = = = = = = = = =
954  ! STARTS LOOP 2 OVER THE SPECTRUM FREQUENCIES
955  ! = = = = = = = = = = = = = = = = = = = = = = = = =
956  DO jf=jfmin,jfmax
957  IF (satval(jf).GT.gqthrsat) THEN
958  !
959  !.........Recovers the coefficient for the coupling factor
960  !.........Computes the coupling coefficients for the case +Delta1 (SIG=1)
961  scal_t=tb_sca(lbuf+jf)*factor
962  t_2p3m=t2p3m*scal_t
963  t_2m3p=t2m3p*scal_t
964  !
965  !.........Frequency indexes and coefficients
966  jfm0=f_poin(jf+lbuf)
967  cf0 =f_coef(jf+lbuf)
968  cp0 =f_proj(jf+lbuf)
969  jfm1=f_poin(jf+ixf1)
970  cf1 =f_coef(jf+ixf1)
971  cp1 =f_proj(jf+ixf1)
972  jfm2=f_poin(jf+ixf2)
973  cf2 =f_coef(jf+ixf2)
974  cp2 =f_proj(jf+ixf2)
975  jfm3=f_poin(jf+ixf3)
976  cf3 =f_coef(jf+ixf3)
977  cp3 =f_proj(jf+ixf3)
978  !
979  ! -------------------------------------------------
980  ! STARTS LOOP 3 OVER THE SPECTRUM DIRECTIONS
981  ! -------------------------------------------------
982  DO jt=1,nt
983  !
984  !...........Direction indexes
985  ! direct config (+delta1) (sig =1)
986  jt1p =t_poin(jt+kt1p)
987  jt1p2p=t_poin(jt+kt1p2p)
988  jt1p2m=t_poin(jt+kt1p2m)
989  jt1p3p=t_poin(jt+kt1p3p)
990  jt1p3m=t_poin(jt+kt1p3m)
991  ! image config (-delta1)
992  jt1m =t_poin(jt+kt1m)
993  jt1m2p=t_poin(jt+kt1m2p)
994  jt1m2m=t_poin(jt+kt1m2m)
995  jt1m3p=t_poin(jt+kt1m3p)
996  jt1m3m=t_poin(jt+kt1m3m)
997  !
998  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
999  ! STARTS LOOP 4 OVER THE MESH NODES
1000  ! - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
1001  !
1002  sp0=f(jt,jfm0)*cf0
1003  !
1004  ! IF (SP0.GT.FSEUIL) THEN
1005  !
1006  ! Config. +Delta1 (SIG=1)
1007  ! =======================
1008  !...............Computes the spectrum values in 1, 2, 3
1009  sp1p =f(jt1p ,jfm1)*cf1
1010  sp1p2p=f(jt1p2p,jfm2)*cf2
1011  sp1p3m=f(jt1p3m,jfm3)*cf3
1012  sp1p2m=f(jt1p2m,jfm2)*cf2
1013  sp1p3p=f(jt1p3p,jfm3)*cf3
1014  !
1015  !...............Computes auxiliary products and variables
1016  aux01=sp0*v1_4+sp1p
1017  aux02=sp0*sp1p
1018  aux03=sp1p2p*sp1p3m
1019  aux04=sp1p2p*v3_4+sp1p3m*v2_4
1020  aux05=sp1p2m*sp1p3p
1021  aux06=sp1p2m*v3_4+sp1p3p*v2_4
1022  aux07=aux02*v3_4
1023  aux08=aux02*v2_4
1024  !
1025  !...............Computes the components of the transfer term
1026  s_2p3m=aux03*aux01-aux02*aux04
1027  s_2m3p=aux05*aux01-aux02*aux06
1028  q_2p3m=t_2p3m*s_2p3m
1029  q_2m3p=t_2m3p*s_2m3p
1030  aux00 =q_2p3m+q_2m3p
1031  !
1032  !...............Computes the components of the derived terms (dQ/dF)
1033  q2pd0 =t_2p3m*(aux03*v1_4 - sp1p*aux04)*cf0
1034  q2pd1 =t_2p3m*(aux03 - sp0 *aux04)*cf1
1035  q2pd2p=t_2p3m*(aux01*sp1p3m - aux07 )*cf2
1036  q2pd3m=t_2p3m*(aux01*sp1p2p - aux08 )*cf3
1037  q2md0 =t_2m3p*(aux05*v1_4 - sp1p*aux06)*cf0
1038  q2md1 =t_2m3p*(aux03 - sp0 *aux06)*cf1
1039  q2md2m=t_2m3p*(aux01*sp1p3p - aux07 )*cf2
1040  q2md3p=t_2m3p*(aux01*sp1p2m - aux08 )*cf3
1041  aux09=q2pd0+q2md0
1042  aux10=q2pd1+q2md1
1043  !
1044  !...............Sum of Qnl4 term in the table TSTOT
1045  tstot(jt,jfm0 )=tstot(jt,jfm0 )+aux00 *cp0
1046  tstot(jt1p,jfm1 )=tstot(jt1p,jfm1 )+aux00 *cp1
1047  tstot(jt1p2p,jfm2)=tstot(jt1p2p,jfm2)-q_2p3m*cp2
1048  tstot(jt1p2m,jfm2)=tstot(jt1p2m,jfm2)-q_2m3p*cp2
1049  tstot(jt1p3m,jfm3)=tstot(jt1p3m,jfm3)-q_2p3m*cp3
1050  tstot(jt1p3p,jfm3)=tstot(jt1p3p,jfm3)-q_2m3p*cp3
1051  !
1052  !...............Sum of the term dQnl4/dF in the table TSDER
1053  tsder(jt,jfm0)=tsder(jt,jfm0)+aux09 *cp0
1054  tsder(jt1p,jfm1)=tsder(jt1p,jfm1)+aux10 *cp1
1055  tsder(jt1p2p,jfm2)=tsder(jt1p2p,jfm2)-q2pd2p*cp2
1056  tsder(jt1p2m,jfm2)=tsder(jt1p2m,jfm2)-q2md2m*cp2
1057  tsder(jt1p3m,jfm3)=tsder(jt1p3m,jfm3)-q2pd3m*cp3
1058  tsder(jt1p3p,jfm3)=tsder(jt1p3p,jfm3)-q2md3p*cp3
1059 #ifdef W3_TGQM
1060  ! Test output to set up triplet method ...
1061  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt, jfm0,aux00 *cp0, f(jt,jfm0),tstot(jt ,jfm0)
1062  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1p, jfm1,aux00 *cp1, f(jt1p,jfm1),tstot(jt1p,jfm1)
1063  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1p2p,jfm2,-q_2p3m*cp2,f(jt1p2p,jfm2),tstot(jt1p2p,jfm2)
1064  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1p2m,jfm2,-q_2m3p*cp2,f(jt1p2m,jfm2),tstot(jt1p2m,jfm2)
1065  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1p3m,jfm2,-q_2p3m*cp3,f(jt1p3m,jfm3),tstot(jt1p3m,jfm3)
1066  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1p3p,jfm2,-q_2m3p*cp3,f(jt1p3p,jfm3),tstot(jt1p3p,jfm3)
1067  temp=(tb_tpm(iq_om2,jt1,jf1)*(( f(jt1p2p,jfm2)*cf2 *f(jt1p3m,jfm3)*cf3)* &
1068  (f(jt,jfm0 )*cf0*tb_v14(jf1)+f(jt1p ,jfm1)*cf1) &
1069  -sp0*sp1p*(sp1p2p*v3_4+sp1p3m*v2_4))+t_2m3p*(aux05*aux01-aux02*aux06)) *cp0
1070  WRITE(995,'(3I3,3E12.3)') iconf,jf,jt, f(jt,jfm0)
1071  temp=(q_2p3m+q_2m3p) *cp1
1072  WRITE(995,'(5I3,3E12.3)') iconf,jf,jt,jt1p, jfm1,aux00 *cp1, f(jt1p,jfm1),tstot(jt1p,jfm1)
1073  WRITE(995,'(5I3,3E12.3)') iconf,jf,jt,jt1p2p,jfm2,-q_2p3m*cp2,f(jt1p2p,jfm2),tstot(jt1p2p,jfm2)
1074  WRITE(995,'(5I3,3E12.3)') iconf,jf,jt,jt1p2m,jfm2,-q_2m3p*cp2,f(jt1p2m,jfm2),tstot(jt1p2m,jfm2)
1075  WRITE(995,'(5I3,3E12.3)') iconf,jf,jt,jt1p3m,jfm2,-q_2p3m*cp3,f(jt1p3m,jfm3),tstot(jt1p3m,jfm3)
1076  WRITE(995,'(5I3,3E12.3)') iconf,jf,jt,jt1p3p,jfm2,-q_2m3p*cp3,f(jt1p3p,jfm3),tstot(jt1p3p,jfm3)
1077 #endif
1078  !
1079  ! Config. -Delta1 (SIG=-1)
1080  ! ========================
1081  !...............Computes the spectrum values in 1, 2, 3
1082  sp1m =f(jt1m ,jfm1)*cf1
1083  sp1m2p=f(jt1m2p,jfm2)*cf2
1084  sp1m3m=f(jt1m3m,jfm3)*cf3
1085  sp1m2m=f(jt1m2m,jfm2)*cf2
1086  sp1m3p=f(jt1m3p,jfm3)*cf3
1087  !
1088  !...............Computes auxiliary products and variables
1089  aux01=sp0*v1_4+sp1m
1090  aux02=sp0*sp1m
1091  aux03=sp1m2p*sp1m3m
1092  aux04=sp1m2p*v3_4+sp1m3m*v2_4
1093  aux05=sp1m2m*sp1m3p
1094  aux06=sp1m2m*v3_4+sp1m3p*v2_4
1095  aux07=aux02*v3_4
1096  aux08=aux02*v2_4
1097  !
1098  !...............Computes the transfer term components
1099  s_2p3m=aux03*aux01-aux02*aux04
1100  s_2m3p=aux05*aux01-aux02*aux06
1101  q_2p3m=t_2m3p*s_2p3m
1102  q_2m3p=t_2p3m*s_2m3p
1103  aux00 =q_2p3m+q_2m3p ! Same as in +Delta1, can be commented out
1104  !
1105  !...............Computes the derived terms components (dQ/dF)
1106  q2pd0 =t_2p3m*(aux03*v1_4 - sp1m*aux04)*cf0
1107  q2pd1 =t_2p3m*(aux03 - sp0 *aux04)*cf1
1108  q2pd2p=t_2p3m*(aux01*sp1m3m - aux07 )*cf2
1109  q2pd3m=t_2p3m*(aux01*sp1m2p - aux08 )*cf3
1110  q2md0 =t_2m3p*(aux05*v1_4 - sp1m*aux06)*cf0
1111  q2md1 =t_2m3p*(aux03 - sp0 *aux06)*cf1
1112  q2md2m=t_2m3p*(aux01*sp1m3p - aux07 )*cf2
1113  q2md3p=t_2m3p*(aux01*sp1m2m - aux08 )*cf3
1114  aux09=q2pd0+q2md0
1115  aux10=q2pd1+q2md1
1116  !
1117  !...............Sum of Qnl4 term in the table TSTOT
1118  tstot(jt ,jfm0)=tstot(jt ,jfm0)+aux00 *cp0
1119  tstot(jt1m ,jfm1)=tstot(jt1m ,jfm1)+aux00 *cp1
1120  tstot(jt1m2p,jfm2)=tstot(jt1m2p,jfm2)-q_2p3m*cp2
1121  tstot(jt1m2m,jfm2)=tstot(jt1m2m,jfm2)-q_2m3p*cp2
1122  tstot(jt1m3m,jfm3)=tstot(jt1m3m,jfm3)-q_2p3m*cp3
1123  tstot(jt1m3p,jfm3)=tstot(jt1m3p,jfm3)-q_2m3p*cp3
1124  !
1125  !...............Sum of the term dQnl4/dF in the table TSDER
1126  tsder(jt ,jfm0)=tsder(jt ,jfm0)+aux09 *cp0
1127  tsder(jt1m ,jfm1)=tsder(jt1m ,jfm1)+aux10 *cp1
1128  tsder(jt1m2p,jfm2)=tsder(jt1m2p,jfm2)-q2pd2p*cp2
1129  tsder(jt1m2m,jfm2)=tsder(jt1m2m,jfm2)-q2md2m*cp2
1130  tsder(jt1m3m,jfm3)=tsder(jt1m3m,jfm3)-q2pd3m*cp3
1131  tsder(jt1m3p,jfm3)=tsder(jt1m3p,jfm3)-q2md3p*cp3
1132  !
1133 #ifdef W3_TGQM
1134  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt, jfm0,aux00 *cp0, f(jt,jfm0),tstot(jt ,jfm0)
1135  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1m, jfm1,aux00 *cp1, f(jt1m,jfm1),tstot(jt1m,jfm1)
1136  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1m2p,jfm2,-q_2p3m*cp2,f(jt1m2p,jfm2),tstot(jt1m2p,jfm2)
1137  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1m2m,jfm2,-q_2m3p*cp2,f(jt1m2m,jfm2),tstot(jt1m2m,jfm2)
1138  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1m3m,jfm2,-q_2p3m*cp3,f(jt1m3m,jfm3),tstot(jt1m3m,jfm3)
1139  WRITE(994,'(5I3,3E12.3)') iconf,jf,jt,jt1m3p,jfm2,-q_2m3p*cp3,f(jt1m3p,jfm3),tstot(jt1m3p,jfm3)
1140 #endif
1141  !
1142  ! ENDIF ! this was the test on SEUIL
1143  !
1144  ENDDO
1145  ! -------------------------------------------------
1146  ! END OF LOOP 3 OVER THE SPECTRUM DIRECTIONS
1147  ! -------------------------------------------------
1148  !
1149  ENDIF ! End of test on saturation level
1150  ENDDO
1151  ! = = = = = = = = = = = = = = = = = = = = = = = = =
1152  ! END OF LOOP 2 OVER THE SPECTRUM FREQUENCIES
1153  ! = = = = = = = = = = = = = = = = = = = = = = = = =
1154  !
1155  ENDDO
1156  ! ==================================================
1157  ! END OF LOOP 1 OVER THE SELECTED CONFIGURATIONS
1158  ! ==================================================
1159  ! Applying WAM DEPTH SCALING ! to be added later ...
1160  ! CALL q_dscale(F,WN,SIG,DTH,NK,NTH,DEPTH,q_dfac)
1161  q_dfac=1
1162 
1163  ! Amplification inspired by Lavrenov 2001, eq 10.
1164  ampfac=gqamp(4)*min(max(accmax/gqamp(2),1.)**gqamp(1),gqamp(3))
1165  !WRITE(991,*) ACCMAX,q_dfac,AMPFAC,GQAMP(1:3),SATVAL(10),SATVAL(30)
1166 
1167  ! Replacing Double Precision with Simple Real and scaling
1168  tstotn = tstot*q_dfac*ampfac
1169  tsdern = tsder*q_dfac*ampfac
1170 
1171 
1172  ! Converting Snl(theta,f) to Snl(theta,k)/sigma
1173  DO ith = 1,nt
1174  DO ik = 1,nf
1175  tstotn(ith,ik) = tstotn(ith,ik)*cg(ik)/(tpi*sig(ik))
1176  ENDDO
1177  ENDDO
1178  !CLOSE(994)
1179  !STOP

References w3gdatmd::dth, f_coef, f_poin, f_proj, w3gdatmd::fr1, w3gdatmd::gqamp, w3gdatmd::gqthrsat, idconf, k_1m, k_1m2m, k_1m2p, k_1m3m, k_1m3p, k_1p, k_1p2m, k_1p2p, k_1p3m, k_1p3p, k_if1, k_if2, k_if3, nconf, w3gdatmd::nk, w3gdatmd::nth, w3gdatmd::sig, t_poin, tb_fac, tb_sca, tb_tmp, tb_tpm, tb_v14, tb_v24, tb_v34, constants::tpi, and w3gdatmd::xfr.

Referenced by w3exnc(), w3expo(), and w3srcemd::w3srce().

Variable Documentation

◆ f_coef

double precision, dimension(:), allocatable w3snl1md::f_coef

Definition at line 93 of file w3snl1md.F90.

93  DOUBLE PRECISION, ALLOCATABLE :: F_COEF(:) , F_PROJ(:) , TB_SCA(:) , TB_V14(:)

Referenced by insnlgqm(), and w3snlgqm().

◆ f_poin

integer, dimension(:), allocatable w3snl1md::f_poin

Definition at line 91 of file w3snl1md.F90.

91  INTEGER, ALLOCATABLE :: F_POIN(:) , T_POIN(:) , K_IF1(:) , K_1P(:,:) , &
92  K_1M(:,:) , IDCONF(:,:)

Referenced by insnlgqm(), and w3snlgqm().

◆ f_proj

double precision, dimension(:), allocatable w3snl1md::f_proj

Definition at line 93 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ idconf

integer, dimension(:,:), allocatable w3snl1md::idconf

Definition at line 91 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1m

integer, dimension(:,:), allocatable w3snl1md::k_1m

Definition at line 91 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1m2m

integer, dimension(:,:,:), allocatable w3snl1md::k_1m2m

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1m2p

integer, dimension(:,:,:), allocatable w3snl1md::k_1m2p

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1m3m

integer, dimension(:,:,:), allocatable w3snl1md::k_1m3m

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1m3p

integer, dimension(:,:,:), allocatable w3snl1md::k_1m3p

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1p

integer, dimension(:,:), allocatable w3snl1md::k_1p

Definition at line 91 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1p2m

integer, dimension(:,:,:), allocatable w3snl1md::k_1p2m

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1p2p

integer, dimension(:,:,:), allocatable w3snl1md::k_1p2p

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1p3m

integer, dimension(:,:,:), allocatable w3snl1md::k_1p3m

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_1p3p

integer, dimension(:,:,:), allocatable w3snl1md::k_1p3p

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_if1

integer, dimension(:), allocatable w3snl1md::k_if1

Definition at line 91 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ k_if2

integer, dimension (:,:,:), allocatable w3snl1md::k_if2

Definition at line 87 of file w3snl1md.F90.

87  INTEGER, ALLOCATABLE :: K_IF2 (:,:,:) , K_IF3 (:,:,:) , K_1P2P(:,:,:) , &
88  K_1P3M(:,:,:) , K_1P2M(:,:,:) , K_1P3P(:,:,:) , &
89  K_1M2P(:,:,:) , K_1M3M(:,:,:) , K_1M2M(:,:,:) , &
90  K_1M3P(:,:,:)

Referenced by insnlgqm(), and w3snlgqm().

◆ k_if3

integer, dimension (:,:,:), allocatable w3snl1md::k_if3

Definition at line 87 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ nconf

integer w3snl1md::nconf

Definition at line 86 of file w3snl1md.F90.

86  INTEGER :: NCONF

Referenced by insnlgqm(), and w3snlgqm().

◆ t_poin

integer, dimension(:), allocatable w3snl1md::t_poin

Definition at line 91 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_fac

double precision, dimension(:,:,:), allocatable w3snl1md::tb_fac

Definition at line 94 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_sca

double precision, dimension(:), allocatable w3snl1md::tb_sca

Definition at line 93 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_tmp

double precision, dimension(:,:,:), allocatable w3snl1md::tb_tmp

Definition at line 94 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_tpm

double precision, dimension(:,:,:), allocatable w3snl1md::tb_tpm

Definition at line 94 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_v14

double precision, dimension(:), allocatable w3snl1md::tb_v14

Definition at line 93 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_v24

double precision, dimension(:,:,:), allocatable w3snl1md::tb_v24

Definition at line 94 of file w3snl1md.F90.

94  DOUBLE PRECISION, ALLOCATABLE :: TB_V24(:,:,:) , TB_V34(:,:,:) , &
95  TB_TPM(:,:,:) , TB_TMP(:,:,:) , TB_FAC(:,:,:)

Referenced by insnlgqm(), and w3snlgqm().

◆ tb_v34

double precision, dimension(:,:,:), allocatable w3snl1md::tb_v34

Definition at line 94 of file w3snl1md.F90.

Referenced by insnlgqm(), and w3snlgqm().

w3gdatmd::nk
integer, pointer nk
Definition: w3gdatmd.F90:1230
w3gdatmd::kdmn
real, pointer kdmn
Definition: w3gdatmd.F90:1347
w3adatmd::awg2
real, pointer awg2
Definition: w3adatmd.F90:666
w3adatmd::swg6
real, pointer swg6
Definition: w3adatmd.F90:666
w3gdatmd::dth
real, pointer dth
Definition: w3gdatmd.F90:1232
w3adatmd::ic52
integer, dimension(:), pointer ic52
Definition: w3adatmd.F90:658
w3adatmd::awg8
real, pointer awg8
Definition: w3adatmd.F90:666
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
w3adatmd::ip14
integer, dimension(:), pointer ip14
Definition: w3adatmd.F90:658
w3adatmd::im13
integer, dimension(:), pointer im13
Definition: w3adatmd.F90:658
w3adatmd::im21
integer, dimension(:), pointer im21
Definition: w3adatmd.F90:658
w3adatmd::ic11
integer, dimension(:), pointer ic11
Definition: w3adatmd.F90:658
w3adatmd::ip13
integer, dimension(:), pointer ip13
Definition: w3adatmd.F90:658
w3adatmd::ic31
integer, dimension(:), pointer ic31
Definition: w3adatmd.F90:658
w3adatmd::ic21
integer, dimension(:), pointer ic21
Definition: w3adatmd.F90:658
w3adatmd::im11
integer, dimension(:), pointer im11
Definition: w3adatmd.F90:658
w3adatmd::swg5
real, pointer swg5
Definition: w3adatmd.F90:666
w3gdatmd::sig
real, dimension(:), pointer sig
Definition: w3gdatmd.F90:1234
w3adatmd::awg6
real, pointer awg6
Definition: w3adatmd.F90:666
w3adatmd::af11
real, dimension(:), pointer af11
Definition: w3adatmd.F90:670
w3adatmd::ic32
integer, dimension(:), pointer ic32
Definition: w3adatmd.F90:658
w3adatmd::ic82
integer, dimension(:), pointer ic82
Definition: w3adatmd.F90:658
w3gdatmd::gqnf1
integer, pointer gqnf1
Definition: w3gdatmd.F90:1345
w3gdatmd::snls2
real, pointer snls2
Definition: w3gdatmd.F90:1347
w3gdatmd::fachfe
real, pointer fachfe
Definition: w3gdatmd.F90:1232
w3adatmd::im22
integer, dimension(:), pointer im22
Definition: w3adatmd.F90:658
w3gdatmd::gqthrsat
real, pointer gqthrsat
Definition: w3gdatmd.F90:1346
w3gdatmd::nltail
real, pointer nltail
Definition: w3gdatmd.F90:1346
w3adatmd::w3dmnl
subroutine w3dmnl(IMOD, NDSE, NDST, NSP, NSPX)
Initialize an individual data grid at the proper dimensions (DIA).
Definition: w3adatmd.F90:2431
w3adatmd::swg4
real, pointer swg4
Definition: w3adatmd.F90:666
w3adatmd::ip21
integer, dimension(:), pointer ip21
Definition: w3adatmd.F90:658
w3arrymd::outmat
subroutine outmat(NDS, A, MX, NX, NY, MNAME)
Definition: w3arrymd.F90:988
w3adatmd::im12
integer, dimension(:), pointer im12
Definition: w3adatmd.F90:658
w3snl4md::twopi
real twopi
Definition: w3snl4md.F90:142
w3adatmd::swg7
real, pointer swg7
Definition: w3adatmd.F90:666
w3adatmd::nfrhgh
integer, pointer nfrhgh
Definition: w3adatmd.F90:657
w3odatmd::ndse
integer, pointer ndse
Definition: w3odatmd.F90:456
scrip_constants::zero
real(kind=scrip_r8), parameter zero
Definition: scrip_constants.f:50
w3adatmd::awg7
real, pointer awg7
Definition: w3adatmd.F90:666
w3adatmd::im23
integer, dimension(:), pointer im23
Definition: w3adatmd.F90:658
w3adatmd::swg2
real, pointer swg2
Definition: w3adatmd.F90:666
w3adatmd::ic72
integer, dimension(:), pointer ic72
Definition: w3adatmd.F90:658
w3servmd
Definition: w3servmd.F90:3
w3gdatmd::gqnt1
integer, pointer gqnt1
Definition: w3gdatmd.F90:1345
w3gdatmd::lam
real, pointer lam
Definition: w3gdatmd.F90:1347
w3adatmd::swg1
real, pointer swg1
Definition: w3adatmd.F90:666
w3adatmd::nspecx
integer, pointer nspecx
Definition: w3adatmd.F90:657
constants::tpiinv
real, parameter tpiinv
TPIINV Inverse of 2*Pi.
Definition: constants.F90:74
w3gdatmd::gqnq_om2
integer, pointer gqnq_om2
Definition: w3gdatmd.F90:1345
w3gdatmd::nth
integer, pointer nth
Definition: w3gdatmd.F90:1230
w3odatmd
Definition: w3odatmd.F90:3
w3adatmd::ip11
integer, dimension(:), pointer ip11
Definition: w3adatmd.F90:658
m_constants::pi
real pi
circular constant, 3.1415...
Definition: mod_constants.f90:29
w3adatmd::ic22
integer, dimension(:), pointer ic22
Definition: w3adatmd.F90:658
w3adatmd::ic12
integer, dimension(:), pointer ic12
Definition: w3adatmd.F90:658
w3adatmd::ic41
integer, dimension(:), pointer ic41
Definition: w3adatmd.F90:658
w3adatmd::ic62
integer, dimension(:), pointer ic62
Definition: w3adatmd.F90:658
w3gdatmd::snls3
real, pointer snls3
Definition: w3gdatmd.F90:1347
w3gdatmd::fr1
real, pointer fr1
Definition: w3gdatmd.F90:1232
w3gdatmd::kdcon
real, pointer kdcon
Definition: w3gdatmd.F90:1347
w3gdatmd::gqamp
real, dimension(:), pointer gqamp
Definition: w3gdatmd.F90:1346
constants::tpi
real, parameter tpi
TPI 2*Pi.
Definition: constants.F90:72
w3adatmd::awg5
real, pointer awg5
Definition: w3adatmd.F90:666
w3adatmd::swg8
real, pointer swg8
Definition: w3adatmd.F90:666
w3servmd::strace
subroutine strace(IENT, SNAME)
Definition: w3servmd.F90:148
w3adatmd::ic71
integer, dimension(:), pointer ic71
Definition: w3adatmd.F90:658
w3adatmd::nfr
integer, pointer nfr
Definition: w3adatmd.F90:657
w3adatmd::ip22
integer, dimension(:), pointer ip22
Definition: w3adatmd.F90:658
w3arrymd
Definition: w3arrymd.F90:3
w3gdatmd::xfr
real, pointer xfr
Definition: w3gdatmd.F90:1232
w3adatmd::awg3
real, pointer awg3
Definition: w3adatmd.F90:666
w3adatmd::ic51
integer, dimension(:), pointer ic51
Definition: w3adatmd.F90:658
w3adatmd::ip12
integer, dimension(:), pointer ip12
Definition: w3adatmd.F90:658
w3gdatmd::snls1
real, pointer snls1
Definition: w3gdatmd.F90:1347
w3adatmd::nspecy
integer, pointer nspecy
Definition: w3adatmd.F90:657
w3odatmd::ndst
integer, pointer ndst
Definition: w3odatmd.F90:456
w3adatmd::im14
integer, dimension(:), pointer im14
Definition: w3adatmd.F90:658
w3adatmd::ic42
integer, dimension(:), pointer ic42
Definition: w3adatmd.F90:658
constants
Define some much-used constants for global use (all defined as PARAMETER).
Definition: constants.F90:20
w3gdatmd
Definition: w3gdatmd.F90:16
w3adatmd::dal2
real, pointer dal2
Definition: w3adatmd.F90:666
w3adatmd::dal1
real, pointer dal1
Definition: w3adatmd.F90:666
w3adatmd::im24
integer, dimension(:), pointer im24
Definition: w3adatmd.F90:658
w3adatmd::ip24
integer, dimension(:), pointer ip24
Definition: w3adatmd.F90:658
w3adatmd::dal3
real, pointer dal3
Definition: w3adatmd.F90:666
w3arrymd::prt2ds
subroutine prt2ds(NDS, NFR0, NFR, NTH, E, FR, UFR, FACSP, FSC, RRCUT, PRVAR, PRUNIT, PNTNME)
Definition: w3arrymd.F90:1943
w3adatmd::nfrchg
integer, pointer nfrchg
Definition: w3adatmd.F90:657
w3adatmd::ic81
integer, dimension(:), pointer ic81
Definition: w3adatmd.F90:658
w3adatmd::awg1
real, pointer awg1
Definition: w3adatmd.F90:666
w3adatmd::swg3
real, pointer swg3
Definition: w3adatmd.F90:666
w3gdatmd::snlc1
real, pointer snlc1
Definition: w3gdatmd.F90:1347
constants::grav
real, parameter grav
GRAV Acc.
Definition: constants.F90:61
w3adatmd::awg4
real, pointer awg4
Definition: w3adatmd.F90:666
w3adatmd::ip23
integer, dimension(:), pointer ip23
Definition: w3adatmd.F90:658
w3gdatmd::gqthrcou
real, pointer gqthrcou
Definition: w3gdatmd.F90:1346
w3adatmd::ic61
integer, dimension(:), pointer ic61
Definition: w3adatmd.F90:658