ROMS
Loading...
Searching...
No Matches
congrad_mod::rprecond Interface Reference

Public Member Functions

subroutine rprecond_dp (ng, lscale, ltrans, outloop, ninnloop, py)
 
subroutine rprecond_r8 (ng, lscale, ltrans, outloop, ninnloop, py)
 

Detailed Description

Definition at line 136 of file congrad.F.

Member Function/Subroutine Documentation

◆ rprecond_dp()

subroutine congrad_mod::rprecond::rprecond_dp ( integer, intent(in) ng,
integer, intent(in) lscale,
logical, intent(in) ltrans,
integer, intent(in) outloop,
integer, intent(in) ninnloop,
real(dp), dimension(ndatum(ng)), intent(inout) py )

Definition at line 1350 of file congrad.F.

1351!***********************************************************************
1352! !
1353! This routine is the preconditioner in observation space. !
1354! !
1355!***********************************************************************
1356!
1357! Imported variable declarations
1358!
1359 logical, intent(in) :: Ltrans
1360!
1361 integer, intent(in) :: ng, Lscale, outLoop, NinnLoop
1362!
1363 real(dp), intent(inout) :: py(Ndatum(ng))
1364!
1365! Local variable declarations.
1366!
1367 integer :: is, ie, inc, iss, i
1368 integer :: nol, nols, nole, ninc
1369 integer :: ingood
1370 integer :: iobs, nvec
1371!
1372 real(dp) :: dla, fac, facritz
1373!
1374! Set the do-loop indices for the sequential preconditioner
1375! loop.
1376!
1377 IF (ltrans) THEN
1378 nols=outloop-1
1379 nole=1
1380 ninc=-1
1381 ELSE
1382 nols=1
1383 nole=outloop-1
1384 ninc=1
1385 END IF
1386!
1387! Apply the preconditioners derived from all previous outer-loops
1388! sequentially.
1389!
1390 DO nol=nols,nole,ninc
1391!
1392! Determine the number of Ritz vectors to use.
1393! For Lritz=.TRUE., choose HvecErr to be larger.
1394!
1395 ingood=0
1396 DO i=1,ninnloop
1397 IF (cg_ritzerr(i,nol).le.ritzmaxerr) THEN
1398 ingood=ingood+1
1399 END IF
1400 END DO
1401 IF (nritzev.gt.0) THEN
1402 nconvritz=nritzev
1403 ingood=nconvritz
1404 ELSE
1405 nconvritz=ingood
1406 END IF
1407!
1408 IF (lscale.gt.0) THEN
1409 is=1
1410 ie=ingood
1411 inc=1
1412 ELSE
1413 is=ingood
1414 ie=1
1415 inc=-1
1416 END IF
1417!
1418 IF (ltrans) THEN
1419 iss=is
1420 is=ie
1421 ie=iss
1422 inc=-inc
1423 END IF
1424!
1425! Lscale determines the form of the preconditioner:
1426!
1427! 1= LMP (NOT CODED YET)
1428! -1= Inverse LMP (NOT CODED YET)
1429! 2= Square root LMP
1430! -2= Inverse square root LMP
1431!
1432 DO nvec=is,ie,inc
1433!
1434 fac=0.0_dp
1435!
1436 IF (lritz) THEN
1437!
1438! Note that the primitive Ritz vectors cg_zv are stored in order of
1439! ascending Ritz values.
1440!
1441 facritz=cg_beta(ninnloop+1,nol)* &
1442 & cg_zv(ninnloop,ninnloop+1-nvec,nol)
1443!
1444! Compute dot-product of py with q_k+1 from previous outer-loop.
1445!
1446 IF (.not.ltrans) THEN
1447 dla=0.0_dp
1448 DO iobs=1,ndatum(ng)
1449 dla=dla+zcglwk(iobs,ninnloop+1,nol)*py(iobs)
1450 END DO
1451 facritz=facritz*dla
1452 END IF
1453
1454 END IF
1455!
1456! Compute the dot-product of py with the Ritz vectors from previous
1457! outer-loop.
1458!
1459! Note that vcglev are arranged in order of descending Ritz values,
1460! while the Ritz values themselves that are stored in cg_Ritz are in
1461! ascending order.
1462!
1463 dla=0.0_dp
1464 DO iobs=1,ndatum(ng)
1465 dla=dla+vcglev(iobs,nvec,nol)*py(iobs)
1466 END DO
1467!
1468! NOTE: Lscale=1 and Lscale=-1 cases are not complete yet.
1469!
1470 IF (lscale.eq.-1) THEN
1471 fac=(cg_ritz(ninnloop+1-nvec,nol)-1.0_dp)*dla
1472 ELSE IF (lscale.eq.1) THEN
1473 fac=(1.0_dp/cg_ritz(ninnloop+1-nvec,nol)-1.0_dp)*dla
1474 ELSE IF (lscale.eq.-2) THEN
1475 fac=(sqrt(cg_ritz(ninnloop+1-nvec,nol))-1.0_dp)*dla
1476 ELSE IF (lscale.eq.2) THEN
1477 fac=(1.0_dp/sqrt(cg_ritz(ninnloop+1-nvec,nol))-1.0_dp)*dla
1478 END IF
1479!
1480 IF (.not.ltrans) THEN
1481 IF (lritz.and.(lscale.eq.-2)) THEN
1482 fac=fac+facritz/sqrt(cg_ritz(ninnloop+1-nvec,nol))
1483 END IF
1484 IF (lritz.and.(lscale.eq.2)) THEN
1485 fac=fac-facritz/cg_ritz(ninnloop+1-nvec,nol)
1486 END IF
1487 END IF
1488!
1489 DO iobs=1,ndatum(ng)
1490 py(iobs)=py(iobs)+fac*vcglev(iobs,nvec,nol)
1491 END DO
1492!
1493 IF (lritz.and.ltrans) THEN
1494
1495 IF (lscale.eq.2) THEN
1496 fac=-facritz*dla/cg_ritz(ninnloop+1-nvec,nol)
1497 END IF
1498 IF (lscale.eq.-2) THEN
1499 fac=facritz*dla/sqrt(cg_ritz(ninnloop+1-nvec,nol))
1500 END IF
1501
1502 DO iobs=1,ndatum(ng)
1503 py(iobs)=py(iobs)+fac*zcglwk(iobs,ninnloop+1,nol)
1504 END DO
1505
1506 END IF
1507
1508 END DO
1509 END DO
1510!
1511 RETURN

References mod_fourdvar::cg_beta, mod_fourdvar::cg_ritz, mod_fourdvar::cg_ritzerr, mod_fourdvar::cg_zv, mod_fourdvar::lritz, mod_fourdvar::nconvritz, mod_fourdvar::nritzev, mod_fourdvar::ritzmaxerr, mod_fourdvar::vcglev, and mod_fourdvar::zcglwk.

◆ rprecond_r8()

subroutine congrad_mod::rprecond::rprecond_r8 ( integer, intent(in) ng,
integer, intent(in) lscale,
logical, intent(in) ltrans,
integer, intent(in) outloop,
integer, intent(in) ninnloop,
real(r8), dimension(ndatum(ng)), intent(inout) py )

Definition at line 1516 of file congrad.F.

1517!***********************************************************************
1518! !
1519! This routine is the preconditioner in observation space. !
1520! !
1521!***********************************************************************
1522!
1523! Imported variable declarations
1524!
1525 logical, intent(in) :: Ltrans
1526!
1527 integer, intent(in) :: ng, Lscale, outLoop, NinnLoop
1528!
1529 real(r8), intent(inout) :: py(Ndatum(ng))
1530!
1531! Local variable declarations.
1532!
1533 integer :: is, ie, inc, iss, i
1534 integer :: nol, nols, nole, ninc
1535 integer :: ingood
1536 integer :: iobs, nvec
1537!
1538 real(r8) :: dla, fac, facritz
1539!
1540! Set the do-loop indices for the sequential preconditioner
1541! loop.
1542!
1543 IF (ltrans) THEN
1544 nols=outloop-1
1545 nole=1
1546 ninc=-1
1547 ELSE
1548 nols=1
1549 nole=outloop-1
1550 ninc=1
1551 END IF
1552!
1553! Apply the preconditioners derived from all previous outer-loops
1554! sequentially.
1555!
1556 DO nol=nols,nole,ninc
1557!
1558! Determine the number of Ritz vectors to use.
1559! For Lritz=.TRUE., choose HvecErr to be larger.
1560!
1561 ingood=0
1562 DO i=1,ninnloop
1563 IF (cg_ritzerr(i,nol).le.ritzmaxerr) THEN
1564 ingood=ingood+1
1565 END IF
1566 END DO
1567 IF (nritzev.gt.0) THEN
1568 nconvritz=nritzev
1569 ingood=nconvritz
1570 ELSE
1571 nconvritz=ingood
1572 END IF
1573!
1574 IF (lscale.gt.0) THEN
1575 is=1
1576 ie=ingood
1577 inc=1
1578 ELSE
1579 is=ingood
1580 ie=1
1581 inc=-1
1582 END IF
1583!
1584 IF (ltrans) THEN
1585 iss=is
1586 is=ie
1587 ie=iss
1588 inc=-inc
1589 END IF
1590!
1591! Lscale determines the form of the preconditioner:
1592!
1593! 1= LMP (NOT CODED YET)
1594! -1= Inverse LMP (NOT CODED YET)
1595! 2= Square root LMP
1596! -2= Inverse square root LMP
1597!
1598 DO nvec=is,ie,inc
1599!
1600 fac=0.0_dp
1601!
1602 IF (lritz) THEN
1603!
1604! Note that the primitive Ritz vectors cg_zv are stored in order of
1605! ascending Ritz values.
1606!
1607 facritz=cg_beta(ninnloop+1,nol)* &
1608 & cg_zv(ninnloop,ninnloop+1-nvec,nol)
1609!
1610! Compute dot-product of py with q_k+1 from previous outer-loop.
1611!
1612 IF (.not.ltrans) THEN
1613 dla=0.0_dp
1614 DO iobs=1,ndatum(ng)
1615 dla=dla+zcglwk(iobs,ninnloop+1,nol)*py(iobs)
1616 END DO
1617 facritz=facritz*dla
1618 END IF
1619
1620 END IF
1621!
1622! Compute the dot-product of py with the Ritz vectors from previous
1623! outer-loop.
1624!
1625! Note that vcglev are arranged in order of descending Ritz values,
1626! while the Ritz values themselves that are stored in cg_Ritz are in
1627! ascending order.
1628!
1629 dla=0.0_dp
1630 DO iobs=1,ndatum(ng)
1631 dla=dla+vcglev(iobs,nvec,nol)*py(iobs)
1632 END DO
1633!
1634! NOTE: Lscale=1 and Lscale=-1 cases are not complete yet.
1635!
1636 IF (lscale.eq.-1) THEN
1637 fac=(cg_ritz(ninnloop+1-nvec,nol)-1.0_dp)*dla
1638 ELSE IF (lscale.eq.1) THEN
1639 fac=(1.0_dp/cg_ritz(ninnloop+1-nvec,nol)-1.0_dp)*dla
1640 ELSE IF (lscale.eq.-2) THEN
1641 fac=(sqrt(cg_ritz(ninnloop+1-nvec,nol))-1.0_dp)*dla
1642 ELSE IF (lscale.eq.2) THEN
1643 fac=(1.0_dp/sqrt(cg_ritz(ninnloop+1-nvec,nol))-1.0_dp)*dla
1644 END IF
1645!
1646 IF (.not.ltrans) THEN
1647 IF (lritz.and.(lscale.eq.-2)) THEN
1648 fac=fac+facritz/sqrt(cg_ritz(ninnloop+1-nvec,nol))
1649 END IF
1650 IF (lritz.and.(lscale.eq.2)) THEN
1651 fac=fac-facritz/cg_ritz(ninnloop+1-nvec,nol)
1652 END IF
1653 END IF
1654!
1655 DO iobs=1,ndatum(ng)
1656 py(iobs)=py(iobs)+fac*vcglev(iobs,nvec,nol)
1657 END DO
1658!
1659 IF (lritz.and.ltrans) THEN
1660
1661 IF (lscale.eq.2) THEN
1662 fac=-facritz*dla/cg_ritz(ninnloop+1-nvec,nol)
1663 END IF
1664 IF (lscale.eq.-2) THEN
1665 fac=facritz*dla/sqrt(cg_ritz(ninnloop+1-nvec,nol))
1666 END IF
1667
1668 DO iobs=1,ndatum(ng)
1669 py(iobs)=py(iobs)+fac*zcglwk(iobs,ninnloop+1,nol)
1670 END DO
1671
1672 END IF
1673
1674 END DO
1675 END DO
1676!
1677 RETURN

References mod_fourdvar::cg_beta, mod_fourdvar::cg_ritz, mod_fourdvar::cg_ritzerr, mod_fourdvar::cg_zv, mod_fourdvar::lritz, mod_fourdvar::nconvritz, mod_fourdvar::nritzev, mod_fourdvar::ritzmaxerr, mod_fourdvar::vcglev, and mod_fourdvar::zcglwk.


The documentation for this interface was generated from the following file: