Skip to content
Snippets Groups Projects
icedyn_rhg_eap.F90 93.3 KiB
Newer Older
Guillaume Samson's avatar
Guillaume Samson committed
1001 1002 1003 1004 1005 1006 1007 1008 1009 1010 1011 1012 1013 1014 1015 1016 1017 1018 1019 1020 1021 1022 1023 1024 1025 1026 1027 1028 1029 1030 1031 1032 1033 1034 1035 1036 1037 1038 1039 1040 1041 1042 1043 1044 1045 1046 1047 1048 1049 1050 1051 1052 1053 1054 1055 1056 1057 1058 1059 1060 1061 1062 1063 1064 1065 1066 1067 1068 1069 1070 1071 1072 1073 1074 1075 1076 1077 1078 1079 1080 1081 1082 1083 1084 1085 1086 1087 1088 1089 1090 1091 1092 1093 1094 1095 1096 1097 1098 1099 1100 1101 1102 1103 1104 1105 1106 1107 1108 1109 1110 1111 1112 1113 1114 1115 1116 1117 1118 1119 1120 1121 1122 1123 1124 1125 1126 1127 1128 1129 1130 1131 1132 1133 1134 1135 1136 1137 1138 1139 1140 1141 1142 1143 1144 1145 1146 1147 1148 1149 1150 1151 1152 1153 1154 1155 1156 1157 1158 1159 1160 1161 1162 1163 1164 1165 1166 1167 1168 1169 1170 1171 1172 1173 1174 1175 1176 1177 1178 1179 1180 1181 1182 1183 1184 1185 1186 1187 1188 1189 1190 1191 1192 1193 1194 1195 1196 1197 1198 1199 1200 1201 1202 1203 1204 1205 1206 1207 1208 1209 1210 1211 1212 1213 1214 1215 1216 1217 1218 1219 1220 1221 1222 1223 1224 1225 1226 1227 1228 1229 1230 1231 1232 1233 1234 1235 1236 1237 1238 1239 1240 1241 1242 1243 1244 1245 1246 1247 1248 1249 1250 1251 1252 1253 1254 1255 1256 1257 1258 1259 1260 1261 1262 1263 1264 1265 1266 1267 1268 1269 1270 1271 1272 1273 1274 1275 1276 1277 1278 1279 1280 1281 1282 1283 1284 1285 1286 1287 1288 1289 1290 1291 1292 1293 1294 1295 1296 1297 1298 1299 1300 1301 1302 1303 1304 1305 1306 1307 1308 1309 1310 1311 1312 1313 1314 1315 1316 1317 1318 1319 1320 1321 1322 1323 1324 1325 1326 1327 1328 1329 1330 1331 1332 1333 1334 1335 1336 1337 1338 1339 1340 1341 1342 1343 1344 1345 1346 1347 1348 1349 1350 1351 1352 1353 1354 1355 1356 1357 1358 1359 1360 1361 1362 1363 1364 1365 1366 1367 1368 1369 1370 1371 1372 1373 1374 1375 1376 1377 1378 1379 1380 1381 1382 1383 1384 1385 1386 1387 1388 1389 1390 1391 1392 1393 1394 1395 1396 1397 1398 1399 1400 1401 1402 1403 1404 1405 1406 1407 1408 1409 1410 1411 1412 1413 1414 1415 1416 1417 1418 1419 1420 1421 1422 1423 1424 1425 1426 1427 1428 1429 1430 1431 1432 1433 1434 1435 1436 1437 1438 1439 1440 1441 1442 1443 1444 1445 1446 1447 1448 1449 1450 1451 1452 1453 1454 1455 1456 1457 1458 1459 1460 1461 1462 1463 1464 1465 1466 1467 1468 1469 1470 1471 1472 1473 1474 1475 1476 1477 1478 1479 1480 1481 1482 1483 1484 1485 1486 1487 1488 1489 1490 1491 1492 1493 1494 1495 1496 1497 1498 1499 1500 1501 1502 1503 1504 1505 1506 1507 1508 1509 1510 1511 1512 1513 1514 1515 1516 1517 1518 1519 1520 1521 1522 1523 1524 1525 1526 1527 1528 1529 1530 1531 1532 1533 1534 1535 1536 1537 1538 1539 1540 1541 1542 1543 1544 1545 1546 1547 1548 1549 1550 1551 1552 1553 1554 1555 1556 1557 1558 1559 1560 1561 1562 1563 1564 1565 1566 1567 1568 1569 1570 1571 1572 1573 1574 1575 1576 1577 1578 1579 1580 1581 1582 1583 1584 1585 1586 1587 1588 1589 1590 1591 1592 1593 1594 1595 1596 1597 1598 1599 1600 1601 1602 1603 1604 1605 1606 1607 1608 1609 1610 1611 1612 1613 1614 1615 1616 1617 1618 1619 1620 1621 1622 1623 1624 1625 1626 1627 1628 1629 1630 1631 1632 1633 1634 1635 1636 1637 1638 1639 1640 1641 1642 1643 1644 1645 1646 1647 1648 1649 1650 1651 1652 1653 1654 1655 1656 1657 1658 1659 1660 1661 1662 1663 1664 1665 1666 1667 1668 1669 1670 1671 1672 1673 1674 1675 1676 1677 1678 1679 1680 1681 1682 1683 1684 1685 1686 1687 1688 1689 1690 1691 1692 1693 1694 1695 1696 1697 1698 1699 1700 1701 1702 1703 1704 1705 1706 1707 1708 1709 1710 1711 1712 1713 1714 1715 1716 1717
      IF( kt == nit000 .AND. kiter == 1 ) THEN
         !
         IF( lwp ) THEN
            WRITE(numout,*)
            WRITE(numout,*) 'rhg_cvg_eap : ice rheology convergence control'
            WRITE(numout,*) '~~~~~~~'
         ENDIF
         !
         IF( lwm ) THEN
            clname = 'ice_cvg.nc'
            IF( .NOT. Agrif_Root() )   clname = TRIM(Agrif_CFixed())//"_"//TRIM(clname)
            istatus = NF90_CREATE( TRIM(clname), NF90_CLOBBER, ncvgid )
            istatus = NF90_DEF_DIM( ncvgid, 'time'  , NF90_UNLIMITED, idtime )
            istatus = NF90_DEF_VAR( ncvgid, 'uice_cvg', NF90_DOUBLE   , (/ idtime /), nvarid )
            istatus = NF90_ENDDEF(ncvgid)
         ENDIF
         !
      ENDIF

      ! time
      it = ( kt - nit000 ) * kitermax + kiter

      ! convergence
      IF( kiter == 1 ) THEN ! remove the first iteration for calculations of convergence (always very large)
         zresm = 0._wp
      ELSE
         zresm = 0._wp
         DO_2D( 0, 0, 0, 0 )
            zresm = MAX( zresm, MAX( ABS( pu(ji,jj) - pub(ji,jj) ) * umask(ji,jj,1), &
               &                     ABS( pv(ji,jj) - pvb(ji,jj) ) * vmask(ji,jj,1) ) * pmsk15(ji,jj) )
         END_2D
         CALL mpp_max( 'icedyn_rhg_evp', zresm )   ! max over the global domain
      ENDIF

      IF( lwm ) THEN
         ! write variables
         istatus = NF90_PUT_VAR( ncvgid, nvarid, (/zresm/), (/it/), (/1/) )
         ! close file
         IF( kt == nitend - nn_fsbc + 1 .AND. kiter == kitermax )   istatus = NF90_CLOSE(ncvgid)
      ENDIF

   END SUBROUTINE rhg_cvg_eap


   SUBROUTINE update_stress_rdg( ksub, kndte, pdivu, ptension, pshear, pa11, pa12, &
      &                          pstressp,  pstressm, pstress12, pstrength, palphar, palphas )
      !!---------------------------------------------------------------------
      !!                   ***  SUBROUTINE update_stress_rdg  ***
      !!
      !! ** Purpose :   Updates the stress depending on values of strain rate and structure
      !!                tensor and for the last subcycle step it computes closing and sliding rate
      !!---------------------------------------------------------------------
      INTEGER,  INTENT(in   ) ::   ksub, kndte
      REAL(wp), INTENT(in   ) ::   pstrength
      REAL(wp), INTENT(in   ) ::   pdivu, ptension, pshear
      REAL(wp), INTENT(in   ) ::   pa11, pa12
      REAL(wp), INTENT(  out) ::   pstressp, pstressm, pstress12
      REAL(wp), INTENT(  out) ::   palphar, palphas

      INTEGER  ::   kx ,ky, ka

      REAL(wp) ::   zstemp11r, zstemp12r, zstemp22r
      REAL(wp) ::   zstemp11s, zstemp12s, zstemp22s
      REAL(wp) ::   za22, zQd11Qd11, zQd11Qd12, zQd12Qd12
      REAL(wp) ::   zQ11Q11, zQ11Q12, zQ12Q12
      REAL(wp) ::   zdtemp11, zdtemp12, zdtemp22
      REAL(wp) ::   zrotstemp11r, zrotstemp12r, zrotstemp22r
      REAL(wp) ::   zrotstemp11s, zrotstemp12s, zrotstemp22s
      REAL(wp) ::   zsig11, zsig12, zsig22
      REAL(wp) ::   zsgprm11, zsgprm12, zsgprm22
      REAL(wp) ::   zAngle_denom_gamma,  zAngle_denom_alpha
      REAL(wp) ::   zTany_1, zTany_2
      REAL(wp) ::   zx, zy, zkxw, kyw, kaw
      REAL(wp) ::   zinvdx, zinvdy, zinvda
      REAL(wp) ::   zdtemp1, zdtemp2, zatempprime

      REAL(wp), PARAMETER ::   ppkfriction = 0.45_wp
      ! Factor to maintain the same stress as in EVP (see Section 3)
      ! Can be set to 1 otherwise
!     REAL(wp), PARAMETER ::   ppinvstressconviso = 1._wp/(1._wp+ppkfriction*ppkfriction)
      REAL(wp), PARAMETER ::   ppinvstressconviso = 1._wp

      ! next statement uses pphi set in main module (icedyn_rhg_eap)
      REAL(wp), PARAMETER ::   ppinvsin = 1._wp/sin(2._wp*pphi) * ppinvstressconviso

      ! compute eigenvalues, eigenvectors and angles for structure tensor, strain
      ! rates

      ! 1) structure tensor
      za22 = 1._wp-pa11
      zQ11Q11 = 1._wp
      zQ12Q12 = rsmall
      zQ11Q12 = rsmall

      ! gamma: angle between general coordiantes and principal axis of A
      ! here Tan2gamma = 2 a12 / (a11 - a22)
      IF((ABS(pa11 - za22) > rsmall).OR.(ABS(pa12) > rsmall)) THEN
         zAngle_denom_gamma = 1._wp/sqrt( ( pa11 - za22 )*( pa11 - za22) + &
                              4._wp*pa12*pa12 )

         zQ11Q11 = 0.5_wp + ( pa11 - za22 )*0.5_wp*zAngle_denom_gamma   !Cos^2
         zQ12Q12 = 0.5_wp - ( pa11 - za22 )*0.5_wp*zAngle_denom_gamma   !Sin^2
         zQ11Q12 = pa12*zAngle_denom_gamma                     !CosSin
      ENDIF

      ! rotation Q*atemp*Q^T
      zatempprime = zQ11Q11*pa11 + 2.0_wp*zQ11Q12*pa12 + zQ12Q12*za22

      ! make first principal value the largest
      zatempprime = max(zatempprime, 1.0_wp - zatempprime)

      ! 2) strain rate
      zdtemp11 = 0.5_wp*(pdivu + ptension)
      zdtemp12 = pshear*0.5_wp
      zdtemp22 = 0.5_wp*(pdivu - ptension)

      ! here Tan2alpha = 2 dtemp12 / (dtemp11 - dtemp22)

      zQd11Qd11 = 1.0_wp
      zQd12Qd12 = rsmall
      zQd11Qd12 = rsmall

      IF((ABS( zdtemp11 - zdtemp22) > rsmall).OR. (ABS(zdtemp12) > rsmall)) THEN

         zAngle_denom_alpha = 1.0_wp/sqrt( ( zdtemp11 - zdtemp22 )* &
                             ( zdtemp11 - zdtemp22 ) + 4.0_wp*zdtemp12*zdtemp12)

         zQd11Qd11 = 0.5_wp + ( zdtemp11 - zdtemp22 )*0.5_wp*zAngle_denom_alpha !Cos^2
         zQd12Qd12 = 0.5_wp - ( zdtemp11 - zdtemp22 )*0.5_wp*zAngle_denom_alpha !Sin^2
         zQd11Qd12 = zdtemp12*zAngle_denom_alpha !CosSin
      ENDIF

      zdtemp1 = zQd11Qd11*zdtemp11 + 2.0_wp*zQd11Qd12*zdtemp12 + zQd12Qd12*zdtemp22
      zdtemp2 = zQd12Qd12*zdtemp11 - 2.0_wp*zQd11Qd12*zdtemp12 + zQd11Qd11*zdtemp22
      ! In cos and sin values
      zx = 0._wp
      IF ((ABS(zdtemp1) > rsmall).OR.(ABS(zdtemp2) > rsmall)) THEN
         zx = atan2(zdtemp2,zdtemp1)
      ENDIF

      ! to ensure the angle lies between pi/4 and 9 pi/4
      IF (zx < rpi*0.25_wp) zx = zx + rpi*2.0_wp

      ! y: angle between major principal axis of strain rate and structure
      ! tensor
      ! y = gamma - alpha
      ! Expressesed componently with
      ! Tany = (Singamma*Cosgamma - Sinalpha*Cosgamma)/(Cos^2gamma - Sin^alpha)

      zTany_1 = zQ11Q12 - zQd11Qd12
      zTany_2 = zQ11Q11 - zQd12Qd12

      zy = 0._wp

      IF ((ABS(zTany_1) > rsmall).OR.(ABS(zTany_2) > rsmall)) THEN
         zy = atan2(zTany_1,zTany_2)
      ENDIF

      ! to make sure y is between 0 and pi
      IF (zy > rpi) zy = zy - rpi
      IF (zy < 0)  zy = zy + rpi

      ! 3) update anisotropic stress tensor
      zinvdx = real(nx_yield-1,kind=wp)/rpi
      zinvdy = real(ny_yield-1,kind=wp)/rpi
      zinvda = 2._wp*real(na_yield-1,kind=wp)

      ! % need 8 coords and 8 weights
      ! % range in kx
      kx  = int((zx-rpi*0.25_wp-rpi)*zinvdx) + 1
      !!clem kx  = MAX( 1, MIN( nx_yield-1, INT((zx-rpi*0.25_wp-rpi)*zinvdx) + 1  ) )
      zkxw = kx - (zx-rpi*0.25_wp-rpi)*zinvdx

      ky  = int(zy*zinvdy) + 1
      !!clem ky  = MAX( 1, MIN( ny_yield-1, INT(zy*zinvdy) + 1 ) )
      kyw = ky - zy*zinvdy

      ka  = int((zatempprime-0.5_wp)*zinvda) + 1
      !!clem ka  = MAX( 1, MIN( na_yield-1, INT((zatempprime-0.5_wp)*zinvda) + 1 ) )
      kaw = ka - (zatempprime-0.5_wp)*zinvda

      ! % Determine sigma_r(A1,Zeta,y) and sigma_s (see Section A1 of Tsamados 2013)
!!$      zstemp11r =  zkxw  * kyw         * kaw         * s11r(kx  ,ky  ,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * kaw         * s11r(kx+1,ky  ,ka  ) &
!!$        & + zkxw         * (1._wp-kyw) * kaw         * s11r(kx  ,ky+1,ka  ) &
!!$        & + zkxw         * kyw         * (1._wp-kaw) * s11r(kx  ,ky  ,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s11r(kx+1,ky+1,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s11r(kx+1,ky  ,ka+1) &
!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s11r(kx  ,ky+1,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s11r(kx+1,ky+1,ka+1)
!!$      zstemp12r =  zkxw  * kyw         * kaw         * s12r(kx  ,ky  ,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * kaw         * s12r(kx+1,ky  ,ka  ) &
!!$        & + zkxw         * (1._wp-kyw) * kaw         * s12r(kx  ,ky+1,ka  ) &
!!$        & + zkxw         * kyw         * (1._wp-kaw) * s12r(kx  ,ky  ,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s12r(kx+1,ky+1,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s12r(kx+1,ky  ,ka+1) &
!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s12r(kx  ,ky+1,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s12r(kx+1,ky+1,ka+1)
!!$      zstemp22r =  zkxw  * kyw         * kaw         * s22r(kx  ,ky  ,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * kaw         * s22r(kx+1,ky  ,ka  ) &
!!$        & + zkxw         * (1._wp-kyw) * kaw         * s22r(kx  ,ky+1,ka  ) &
!!$        & + zkxw         * kyw         * (1._wp-kaw) * s22r(kx  ,ky  ,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s22r(kx+1,ky+1,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s22r(kx+1,ky  ,ka+1) &
!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s22r(kx  ,ky+1,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s22r(kx+1,ky+1,ka+1)
!!$
!!$      zstemp11s =  zkxw  * kyw         * kaw         * s11s(kx  ,ky  ,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * kaw         * s11s(kx+1,ky  ,ka  ) &
!!$        & + zkxw         * (1._wp-kyw) * kaw         * s11s(kx  ,ky+1,ka  ) &
!!$        & + zkxw         * kyw         * (1._wp-kaw) * s11s(kx  ,ky  ,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s11s(kx+1,ky+1,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s11s(kx+1,ky  ,ka+1) &
!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s11s(kx  ,ky+1,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s11s(kx+1,ky+1,ka+1)
!!$      zstemp12s =  zkxw  * kyw         * kaw         * s12s(kx  ,ky  ,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * kaw         * s12s(kx+1,ky  ,ka  ) &
!!$        & + zkxw         * (1._wp-kyw) * kaw         * s12s(kx  ,ky+1,ka  ) &
!!$        & + zkxw         * kyw         * (1._wp-kaw) * s12s(kx  ,ky  ,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s12s(kx+1,ky+1,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s12s(kx+1,ky  ,ka+1) &
!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s12s(kx  ,ky+1,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s12s(kx+1,ky+1,ka+1)
!!$      zstemp22s =  zkxw  * kyw         * kaw         * s22s(kx  ,ky  ,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * kaw         * s22s(kx+1,ky  ,ka  ) &
!!$        & + zkxw         * (1._wp-kyw) * kaw         * s22s(kx  ,ky+1,ka  ) &
!!$        & + zkxw         * kyw         * (1._wp-kaw) * s22s(kx  ,ky  ,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * kaw         * s22s(kx+1,ky+1,ka  ) &
!!$        & + (1._wp-zkxw) * kyw         * (1._wp-kaw) * s22s(kx+1,ky  ,ka+1) &
!!$        & + zkxw         * (1._wp-kyw) * (1._wp-kaw) * s22s(kx  ,ky+1,ka+1) &
!!$        & + (1._wp-zkxw) * (1._wp-kyw) * (1._wp-kaw) * s22s(kx+1,ky+1,ka+1)

      zstemp11r = s11r(kx,ky,ka)
      zstemp12r = s12r(kx,ky,ka)
      zstemp22r = s22r(kx,ky,ka)
      zstemp11s = s11s(kx,ky,ka)
      zstemp12s = s12s(kx,ky,ka)
      zstemp22s = s22s(kx,ky,ka)


      ! Calculate mean ice stress over a collection of floes (Equation 3 in
      ! Tsamados 2013)

      zsig11  = pstrength*(zstemp11r + ppkfriction*zstemp11s) * ppinvsin
      zsig12  = pstrength*(zstemp12r + ppkfriction*zstemp12s) * ppinvsin
      zsig22  = pstrength*(zstemp22r + ppkfriction*zstemp22s) * ppinvsin

      ! Back - rotation of the stress from principal axes into general coordinates

      ! Update stress
      zsgprm11 = zQ11Q11*zsig11 + zQ12Q12*zsig22 -      2._wp*zQ11Q12 *zsig12
      zsgprm12 = zQ11Q12*zsig11 - zQ11Q12*zsig22 + (zQ11Q11 - zQ12Q12)*zsig12
      zsgprm22 = zQ12Q12*zsig11 + zQ11Q11*zsig22 +      2._wp*zQ11Q12 *zsig12

      pstressp  = zsgprm11 + zsgprm22
      pstress12 = zsgprm12
      pstressm  = zsgprm11 - zsgprm22

      ! Compute ridging and sliding functions in general coordinates
      ! (Equation 11 in Tsamados 2013)
      IF (ksub == kndte) THEN
         zrotstemp11r = zQ11Q11*zstemp11r - 2._wp*zQ11Q12* zstemp12r &
                      + zQ12Q12*zstemp22r
         zrotstemp12r = zQ11Q11*zstemp12r +       zQ11Q12*(zstemp11r-zstemp22r) &
                      - zQ12Q12*zstemp12r
         zrotstemp22r = zQ12Q12*zstemp11r + 2._wp*zQ11Q12* zstemp12r &
                      + zQ11Q11*zstemp22r

         zrotstemp11s = zQ11Q11*zstemp11s - 2._wp*zQ11Q12* zstemp12s &
                      + zQ12Q12*zstemp22s
         zrotstemp12s = zQ11Q11*zstemp12s +       zQ11Q12*(zstemp11s-zstemp22s) &
                      - zQ12Q12*zstemp12s
         zrotstemp22s = zQ12Q12*zstemp11s + 2._wp*zQ11Q12* zstemp12s &
                      + zQ11Q11*zstemp22s

         palphar =  zrotstemp11r*zdtemp11 + 2._wp*zrotstemp12r*zdtemp12 &
                  + zrotstemp22r*zdtemp22
         palphas =  zrotstemp11s*zdtemp11 + 2._wp*zrotstemp12s*zdtemp12 &
                  + zrotstemp22s*zdtemp22
      ENDIF
   END SUBROUTINE update_stress_rdg

!=======================================================================


   SUBROUTINE calc_ffrac( pstressp, pstressm, pstress12, pa11, pa12, &
      &                   pmresult11, pmresult12 )
      !!---------------------------------------------------------------------
      !!                   ***  ROUTINE calc_ffrac  ***
      !!
      !! ** Purpose :   Computes term in evolution equation for structure tensor
      !!                which determines the ice floe re-orientation due to fracture
      !! ** Method :    Eq. 7: Ffrac = -kf(A-S) or = 0 depending on sigma_1 and sigma_2
      !!---------------------------------------------------------------------
      REAL (wp), INTENT(in)  ::   pstressp, pstressm, pstress12, pa11, pa12
      REAL (wp), INTENT(out) ::   pmresult11, pmresult12

      ! local variables

      REAL (wp) ::   zsigma11, zsigma12, zsigma22  ! stress tensor elements
      REAL (wp) ::   zAngle_denom        ! angle with principal component axis
      REAL (wp) ::   zsigma_1, zsigma_2   ! principal components of stress
      REAL (wp) ::   zQ11, zQ12, zQ11Q11, zQ11Q12, zQ12Q12

!!$   REAL (wp), PARAMETER ::   ppkfrac = 0.0001_wp   ! rate of fracture formation
      REAL (wp), PARAMETER ::   ppkfrac = 1.e-3_wp    ! rate of fracture formation
      REAL (wp), PARAMETER ::   ppthreshold = 0.3_wp  ! critical confinement ratio
      !!---------------------------------------------------------------
      !
      zsigma11 = 0.5_wp*(pstressp+pstressm)
      zsigma12 = pstress12
      zsigma22 = 0.5_wp*(pstressp-pstressm)

      ! Here's the change - no longer calculate gamma,
      ! use trig formulation, angle signs are kept correct, don't worry

      ! rotate tensor to get into sigma principal axis

      ! here Tan2gamma = 2 sig12 / (sig11 - sig12)
      ! add rsmall to the denominator to stop 1/0 errors, makes very little
      ! error to the calculated angles < rsmall

      zQ11Q11 = 0.1_wp
      zQ12Q12 = rsmall
      zQ11Q12 = rsmall

      IF((ABS( zsigma11 - zsigma22) > rsmall).OR.(ABS(zsigma12) > rsmall)) THEN

         zAngle_denom = 1.0_wp/sqrt( ( zsigma11 - zsigma22 )*( zsigma11 - &
                       zsigma22 ) + 4.0_wp*zsigma12*zsigma12)

         zQ11Q11 = 0.5_wp + ( zsigma11 - zsigma22 )*0.5_wp*zAngle_denom   !Cos^2
         zQ12Q12 = 0.5_wp - ( zsigma11 - zsigma22 )*0.5_wp*zAngle_denom   !Sin^2
         zQ11Q12 = zsigma12*zAngle_denom                      !CosSin
      ENDIF

      zsigma_1 = zQ11Q11*zsigma11 + 2.0_wp*zQ11Q12*zsigma12 + zQ12Q12*zsigma22 ! S(1,1)
      zsigma_2 = zQ12Q12*zsigma11 - 2.0_wp*zQ11Q12*zsigma12 + zQ11Q11*zsigma22 ! S(2,2)

      ! Pure divergence
      IF ((zsigma_1 >= 0.0_wp).AND.(zsigma_2 >= 0.0_wp))  THEN
         pmresult11 = 0.0_wp
         pmresult12 = 0.0_wp

      ! Unconfined compression: cracking of blocks not along the axial splitting
      ! direction
      ! which leads to the loss of their shape, so we again model it through diffusion
      ELSEIF ((zsigma_1 >= 0.0_wp).AND.(zsigma_2 < 0.0_wp))  THEN
         pmresult11 = - ppkfrac * (pa11 - zQ12Q12)
         pmresult12 = - ppkfrac * (pa12 + zQ11Q12)

      ! Shear faulting
      ELSEIF (zsigma_2 == 0.0_wp) THEN
         pmresult11 = 0.0_wp
         pmresult12 = 0.0_wp
      ELSEIF ((zsigma_1 <= 0.0_wp).AND.(zsigma_1/zsigma_2 <= ppthreshold)) THEN
         pmresult11 = - ppkfrac * (pa11 - zQ12Q12)
         pmresult12 = - ppkfrac * (pa12 + zQ11Q12)

      ! Horizontal spalling
      ELSE
         pmresult11 = 0.0_wp
         pmresult12 = 0.0_wp
      ENDIF

   END SUBROUTINE calc_ffrac


   SUBROUTINE rhg_eap_rst( cdrw, kt )
      !!---------------------------------------------------------------------
      !!                   ***  ROUTINE rhg_eap_rst  ***
      !!
      !! ** Purpose :   Read or write RHG file in restart file
      !!
      !! ** Method  :   use of IOM library
      !!----------------------------------------------------------------------
      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
      !
      INTEGER  ::   iter                      ! local integer
      INTEGER  ::   id1, id2, id3, id4, id5   ! local integers
      INTEGER  ::   ix, iy, ip, iz, n, ia     ! local integers

      INTEGER, PARAMETER            ::    nz = 100

      REAL(wp) ::   ainit, xinit, yinit, pinit, zinit
      REAL(wp) ::   da, dx, dy, dp, dz, a1

      !!clem
      REAL(wp) ::   zw1, zw2, zfac, ztemp
      REAL(wp) ::   zidx, zidy, zidz
      REAL(wp) ::   zsaak(6)                  ! temporary array

      REAL(wp), PARAMETER           ::   eps6 = 1.0e-6_wp
      !!----------------------------------------------------------------------
      !
      IF( TRIM(cdrw) == 'READ' ) THEN        ! Read/initialize
         !                                   ! ---------------
         IF( ln_rstart ) THEN                   !* Read the restart file
            !
            id1 = iom_varid( numrir, 'stress1_i' , ldstop = .FALSE. )
            id2 = iom_varid( numrir, 'stress2_i' , ldstop = .FALSE. )
            id3 = iom_varid( numrir, 'stress12_i', ldstop = .FALSE. )
            id4 = iom_varid( numrir, 'aniso_11'  , ldstop = .FALSE. )
            id5 = iom_varid( numrir, 'aniso_12'  , ldstop = .FALSE. )
            !
            IF( MIN( id1, id2, id3, id4, id5 ) > 0 ) THEN      ! fields exist
               CALL iom_get( numrir, jpdom_auto, 'stress1_i' , stress1_i , cd_type = 'T' )
               CALL iom_get( numrir, jpdom_auto, 'stress2_i' , stress2_i , cd_type = 'T' )
               CALL iom_get( numrir, jpdom_auto, 'stress12_i', stress12_i, cd_type = 'F' )
               CALL iom_get( numrir, jpdom_auto, 'aniso_11'  , aniso_11  , cd_type = 'T' )
               CALL iom_get( numrir, jpdom_auto, 'aniso_12'  , aniso_12  , cd_type = 'T' )
            ELSE                                     ! start rheology from rest
               IF(lwp) WRITE(numout,*)
               IF(lwp) WRITE(numout,*) '   ==>>>   previous run without rheology, set stresses to 0'
               stress1_i (:,:) = 0._wp
               stress2_i (:,:) = 0._wp
               stress12_i(:,:) = 0._wp
               aniso_11  (:,:) = 0.5_wp
               aniso_12  (:,:) = 0._wp
            ENDIF
         ELSE                                   !* Start from rest
            IF(lwp) WRITE(numout,*)
            IF(lwp) WRITE(numout,*) '   ==>>>   start from rest: set stresses to 0'
            stress1_i (:,:) = 0._wp
            stress2_i (:,:) = 0._wp
            stress12_i(:,:) = 0._wp
            aniso_11  (:,:) = 0.5_wp
            aniso_12  (:,:) = 0._wp
         ENDIF
         !

         da = 0.5_wp/real(na_yield-1,kind=wp)
         ainit = 0.5_wp - da
         dx = rpi/real(nx_yield-1,kind=wp)
         xinit = rpi + 0.25_wp*rpi - dx
         dz = rpi/real(nz,kind=wp)
         zinit = -rpi*0.5_wp
         dy = rpi/real(ny_yield-1,kind=wp)
         yinit = -dy

         s11r(:,:,:) = 0._wp
         s12r(:,:,:) = 0._wp
         s22r(:,:,:) = 0._wp
         s11s(:,:,:) = 0._wp
         s12s(:,:,:) = 0._wp
         s22s(:,:,:) = 0._wp

!!$         DO ia=1,na_yield
!!$            DO ix=1,nx_yield
!!$               DO iy=1,ny_yield
!!$                  s11r(ix,iy,ia) = 0._wp
!!$                  s12r(ix,iy,ia) = 0._wp
!!$                  s22r(ix,iy,ia) = 0._wp
!!$                  s11s(ix,iy,ia) = 0._wp
!!$                  s12s(ix,iy,ia) = 0._wp
!!$                  s22s(ix,iy,ia) = 0._wp
!!$                  IF (ia <= na_yield-1) THEN
!!$                     DO iz=1,nz
!!$                        s11r(ix,iy,ia) = s11r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
!!$                           s11kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
!!$                        s12r(ix,iy,ia) = s12r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
!!$                           s12kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
!!$                        s22r(ix,iy,ia) = s22r(ix,iy,ia) + 1*w1(ainit+ia*da)* &
!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
!!$                           s22kr(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
!!$                        s11s(ix,iy,ia) = s11s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
!!$                           s11ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
!!$                        s12s(ix,iy,ia) = s12s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
!!$                           s12ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
!!$                        s22s(ix,iy,ia) = s22s(ix,iy,ia) + 1*w1(ainit+ia*da)* &
!!$                           exp(-w2(ainit+ia*da)*(zinit+iz*dz)*(zinit+iz*dz))* &
!!$                           s22ks(xinit+ix*dx,yinit+iy*dy,zinit+iz*dz)*dz/sin(2._wp*pphi)
!!$                     ENDDO
!!$                     IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
!!$                     IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
!!$                     IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
!!$                     IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
!!$                     IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
!!$                     IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
!!$                  ELSE
!!$                     s11r(ix,iy,ia) = 0.5_wp*s11kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
!!$                     s12r(ix,iy,ia) = 0.5_wp*s12kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
!!$                     s22r(ix,iy,ia) = 0.5_wp*s22kr(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
!!$                     s11s(ix,iy,ia) = 0.5_wp*s11ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
!!$                     s12s(ix,iy,ia) = 0.5_wp*s12ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
!!$                     s22s(ix,iy,ia) = 0.5_wp*s22ks(xinit+ix*dx,yinit+iy*dy,0._wp)/sin(2._wp*pphi)
!!$                     IF (abs(s11r(ix,iy,ia)) < eps6) s11r(ix,iy,ia) = 0._wp
!!$                     IF (abs(s12r(ix,iy,ia)) < eps6) s12r(ix,iy,ia) = 0._wp
!!$                     IF (abs(s22r(ix,iy,ia)) < eps6) s22r(ix,iy,ia) = 0._wp
!!$                     IF (abs(s11s(ix,iy,ia)) < eps6) s11s(ix,iy,ia) = 0._wp
!!$                     IF (abs(s12s(ix,iy,ia)) < eps6) s12s(ix,iy,ia) = 0._wp
!!$                     IF (abs(s22s(ix,iy,ia)) < eps6) s22s(ix,iy,ia) = 0._wp
!!$                  ENDIF
!!$               ENDDO
!!$            ENDDO
!!$         ENDDO

         !! faster but still very slow => to be improved
         zfac = dz/sin(2._wp*pphi)
         DO ia = 1, na_yield-1
            zw1 = w1(ainit+ia*da)
            zw2 = w2(ainit+ia*da)
            DO iz = 1, nz
               zidz = zinit+iz*dz
               ztemp = zw1 * EXP(-zw2*(zinit+iz*dz)*(zinit+iz*dz))
               DO iy = 1, ny_yield
                  zidy = yinit+iy*dy
                  DO ix = 1, nx_yield
                     zidx = xinit+ix*dx
                     call all_skr_sks(zidx,zidy,zidz,zsaak)
                     zsaak = ztemp*zsaak*zfac
                     s11r(ix,iy,ia) = s11r(ix,iy,ia) + zsaak(1)
                     s12r(ix,iy,ia) = s12r(ix,iy,ia) + zsaak(2)
                     s22r(ix,iy,ia) = s22r(ix,iy,ia) + zsaak(3)
                     s11s(ix,iy,ia) = s11s(ix,iy,ia) + zsaak(4)
                     s12s(ix,iy,ia) = s12s(ix,iy,ia) + zsaak(5)
                     s22s(ix,iy,ia) = s22s(ix,iy,ia) + zsaak(6)
                  END DO
               END DO
            END DO
         END DO
         zfac = 1._wp/sin(2._wp*pphi)
         ia = na_yield
         DO iy = 1, ny_yield
            zidy = yinit+iy*dy
            DO ix = 1, nx_yield
               zidx = xinit+ix*dx
               call all_skr_sks(zidx,zidy,0._wp,zsaak)
               zsaak = 0.5_wp*zsaak*zfac
               s11r(ix,iy,ia) = zsaak(1)
               s12r(ix,iy,ia) = zsaak(2)
               s22r(ix,iy,ia) = zsaak(3)
               s11s(ix,iy,ia) = zsaak(4)
               s12s(ix,iy,ia) = zsaak(5)
               s22s(ix,iy,ia) = zsaak(6)
            ENDDO
         ENDDO
         WHERE (ABS(s11r(:,:,:)) < eps6) s11r(:,:,:) = 0._wp
         WHERE (ABS(s12r(:,:,:)) < eps6) s12r(:,:,:) = 0._wp
         WHERE (ABS(s22r(:,:,:)) < eps6) s22r(:,:,:) = 0._wp
         WHERE (ABS(s11s(:,:,:)) < eps6) s11s(:,:,:) = 0._wp
         WHERE (ABS(s12s(:,:,:)) < eps6) s12s(:,:,:) = 0._wp
         WHERE (ABS(s22s(:,:,:)) < eps6) s22s(:,:,:) = 0._wp


      ELSEIF( TRIM(cdrw) == 'WRITE' ) THEN   ! Create restart file
         !                                   ! -------------------
         IF(lwp) WRITE(numout,*) '---- rhg-rst ----'
         iter = kt + nn_fsbc - 1             ! ice restarts are written at kt == nitrst - nn_fsbc + 1
         !
         CALL iom_rstput( iter, nitrst, numriw, 'stress1_i' , stress1_i  )
         CALL iom_rstput( iter, nitrst, numriw, 'stress2_i' , stress2_i  )
         CALL iom_rstput( iter, nitrst, numriw, 'stress12_i', stress12_i )
         CALL iom_rstput( iter, nitrst, numriw, 'aniso_11'  , aniso_11   )
         CALL iom_rstput( iter, nitrst, numriw, 'aniso_12'  , aniso_12   )
         !
      ENDIF
      !
   END SUBROUTINE rhg_eap_rst


   FUNCTION w1(pa)
      !!-------------------------------------------------------------------
      !! Function : w1 (see Gaussian function psi in Tsamados et al 2013)
      !!-------------------------------------------------------------------
      REAL(wp), INTENT(in   ) ::   pa
      REAL(wp) ::   w1
      !!-------------------------------------------------------------------

      w1 = -   223.87569446_wp &
       &   +  2361.21986630_wp*pa &
       &   - 10606.56079975_wp*pa*pa &
       &   + 26315.50025642_wp*pa*pa*pa &
       &   - 38948.30444297_wp*pa*pa*pa*pa &
       &   + 34397.72407466_wp*pa*pa*pa*pa*pa &
       &   - 16789.98003081_wp*pa*pa*pa*pa*pa*pa &
       &   +  3495.82839237_wp*pa*pa*pa*pa*pa*pa*pa

   END FUNCTION w1


   FUNCTION w2(pa)
      !!-------------------------------------------------------------------
      !! Function : w2 (see Gaussian function psi in Tsamados et al 2013)
      !!-------------------------------------------------------------------
      REAL(wp), INTENT(in   ) ::   pa
      REAL(wp) ::   w2
      !!-------------------------------------------------------------------

      w2 = -    6670.68911883_wp &
       &   +   70222.33061536_wp*pa &
       &   -  314871.71525448_wp*pa*pa &
       &   +  779570.02793492_wp*pa*pa*pa &
       &   - 1151098.82436864_wp*pa*pa*pa*pa &
       &   + 1013896.59464498_wp*pa*pa*pa*pa*pa &
       &   -  493379.44906738_wp*pa*pa*pa*pa*pa*pa &
       &   +  102356.55151800_wp*pa*pa*pa*pa*pa*pa*pa

   END FUNCTION w2

   SUBROUTINE all_skr_sks( px, py, pz, allsk )
      REAL(wp), INTENT(in   ) ::   px,py,pz
      REAL(wp), INTENT(out  ) ::   allsk(6)

      REAL(wp) ::   zs12r0, zs21r0
      REAL(wp) ::   zs12s0, zs21s0

      REAL(wp) ::   zpih
      REAL(wp) ::   zn1t2i11, zn1t2i12, zn1t2i21, zn1t2i22
      REAL(wp) ::   zn2t1i11, zn2t1i12, zn2t1i21, zn2t1i22
      REAL(wp) ::   zt1t2i11, zt1t2i12, zt1t2i21, zt1t2i22
      REAL(wp) ::   zt2t1i11, zt2t1i12, zt2t1i21, zt2t1i22
      REAL(wp) ::   zd11, zd12, zd22
      REAL(wp) ::   zIIn1t2, zIIn2t1, zIIt1t2
      REAL(wp) ::   zHen1t2, zHen2t1
      !!-------------------------------------------------------------------

      zpih = 0.5_wp*rpi

      zn1t2i11 = cos(pz+zpih-pphi) * cos(pz+pphi)
      zn1t2i12 = cos(pz+zpih-pphi) * sin(pz+pphi)
      zn1t2i21 = sin(pz+zpih-pphi) * cos(pz+pphi)
      zn1t2i22 = sin(pz+zpih-pphi) * sin(pz+pphi)
      zn2t1i11 = cos(pz-zpih+pphi) * cos(pz-pphi)
      zn2t1i12 = cos(pz-zpih+pphi) * sin(pz-pphi)
      zn2t1i21 = sin(pz-zpih+pphi) * cos(pz-pphi)
      zn2t1i22 = sin(pz-zpih+pphi) * sin(pz-pphi)
      zt1t2i11 = cos(pz-pphi) * cos(pz+pphi)
      zt1t2i12 = cos(pz-pphi) * sin(pz+pphi)
      zt1t2i21 = sin(pz-pphi) * cos(pz+pphi)
      zt1t2i22 = sin(pz-pphi) * sin(pz+pphi)
      zt2t1i11 = cos(pz+pphi) * cos(pz-pphi)
      zt2t1i12 = cos(pz+pphi) * sin(pz-pphi)
      zt2t1i21 = sin(pz+pphi) * cos(pz-pphi)
      zt2t1i22 = sin(pz+pphi) * sin(pz-pphi)
   ! In expression of tensor d, with this formulatin d(x)=-d(x+pi)
   ! Solution, when diagonalizing always check sgn(a11-a22) if > then keep x else
   ! x=x-pi/2
      zd11 = cos(py)*cos(py)*(cos(px)+sin(px)*tan(py)*tan(py))
      zd12 = cos(py)*cos(py)*tan(py)*(-cos(px)+sin(px))
      zd22 = cos(py)*cos(py)*(sin(px)+cos(px)*tan(py)*tan(py))
      zIIn1t2 = zn1t2i11 * zd11 + (zn1t2i12 + zn1t2i21) * zd12 + zn1t2i22 * zd22
      zIIn2t1 = zn2t1i11 * zd11 + (zn2t1i12 + zn2t1i21) * zd12 + zn2t1i22 * zd22
      zIIt1t2 = zt1t2i11 * zd11 + (zt1t2i12 + zt1t2i21) * zd12 + zt1t2i22 * zd22

      IF (-zIIn1t2>=rsmall) THEN
      zHen1t2 = 1._wp
      ELSE
      zHen1t2 = 0._wp
      ENDIF

      IF (-zIIn2t1>=rsmall) THEN
      zHen2t1 = 1._wp
      ELSE
      zHen2t1 = 0._wp
      ENDIF

      !!-------------------------------------------------------------------
      !! Function : s11kr
      !!-------------------------------------------------------------------
      allsk(1) = (- zHen1t2 * zn1t2i11 - zHen2t1 * zn2t1i11)
      !!-------------------------------------------------------------------
      !! Function : s12kr
      !!-------------------------------------------------------------------
      zs12r0 = (- zHen1t2 * zn1t2i12 - zHen2t1 * zn2t1i12)
      zs21r0 = (- zHen1t2 * zn1t2i21 - zHen2t1 * zn2t1i21)
      allsk(2)=0.5_wp*(zs12r0+zs21r0)
      !!-------------------------------------------------------------------
      !! Function : s22kr
      !!-------------------------------------------------------------------
      allsk(3) = (- zHen1t2 * zn1t2i22 - zHen2t1 * zn2t1i22)
      !!-------------------------------------------------------------------
      !! Function : s11ks
      !!-------------------------------------------------------------------

      allsk(4) = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i11 + zHen2t1 * zt2t1i11)
      !!-------------------------------------------------------------------
      !! Function : s12ks
      !!-------------------------------------------------------------------
      zs12s0 = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i12 + zHen2t1 * zt2t1i12)
      zs21s0 = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i21 + zHen2t1 * zt2t1i21)
      allsk(5)=0.5_wp*(zs12s0+zs21s0)
      !!-------------------------------------------------------------------
      !! Function : s22ks
      !!-------------------------------------------------------------------
      allsk(6) = sign(1._wp,zIIt1t2+rsmall)*(zHen1t2 * zt1t2i22 + zHen2t1 * zt2t1i22)
   END SUBROUTINE all_skr_sks

#else
   !!----------------------------------------------------------------------
   !!   Default option         Empty module           NO SI3 sea-ice model
   !!----------------------------------------------------------------------
   USE par_kind
   USE lib_mpp
   CONTAINS
   SUBROUTINE ice_dyn_rhg_eap( kt, Kmm, pstress1_i, pstress2_i, pstress12_i, pshear_i, pdivu_i, pdelta_i, paniso_11, paniso_12, prdg_conv )
      INTEGER                 , INTENT(in   ) ::   kt                                    ! time step
      INTEGER                 , INTENT(in   ) ::   Kmm                                   ! ocean time level index
      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   pstress1_i, pstress2_i, pstress12_i   !
      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   pshear_i  , pdivu_i   , pdelta_i      !
      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   paniso_11 , paniso_12                 ! structure tensor components
      REAL(wp), DIMENSION(:,:), INTENT(in   ) ::   prdg_conv                             ! for ridging
      CALL ctl_stop('EAP rheology is currently dsabled due to issues with AGRIF preprocessing')
   END SUBROUTINE ice_dyn_rhg_eap
   SUBROUTINE rhg_eap_rst( cdrw, kt )
      CHARACTER(len=*) , INTENT(in) ::   cdrw   ! "READ"/"WRITE" flag
      INTEGER, OPTIONAL, INTENT(in) ::   kt     ! ice time-step
   END SUBROUTINE rhg_eap_rst
#endif

   !!==============================================================================
END MODULE icedyn_rhg_eap