diva.f 3.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126
  1. ! ==================================================================
  2. ! ------------------------------------------------------------------
  3. !
  4. subroutine diva
  5. use lsgvar
  6. implicit none
  7. !
  8. ! ------------------------------------------------------------------
  9. !
  10. !**** *diva*.
  11. !
  12. ! by E. Maier-Reimer.
  13. ! Last modified by U. Mikolajewicz 9/87.
  14. !
  15. ! Purpose.
  16. ! --------
  17. ! *diva* controls the computation of diagnostic variables
  18. ! and computes total velocities (*utot* and *vtot*).
  19. ! Dependent on *nsve* different subroutines are called.
  20. !
  21. !
  22. ! *nsve*=2 (exact solution for barotropic velocities).
  23. ! --------
  24. ! *diva* calls *press* computes normalized pressure *p*.
  25. ! *uvtrop2* computes barotropic velocities *ub*,
  26. ! *vb* and surface elevation *zeta*.
  27. ! *geost2* computes baroclinic velocities.
  28. ! *cont1* computes vertical velocities *w*.
  29. !
  30. !** Input.
  31. ! ------
  32. ! rhdif pot density-difference in the vertical direction.
  33. ! taux ) windstress divided by density.
  34. ! tauy )
  35. ! old velocities and surface elevation.
  36. ! all via common /lsgfie/.
  37. !
  38. ! Output.
  39. ! -------
  40. ! utot ) total velocities.
  41. ! vtot )
  42. ! w vertical velocities.
  43. ! ub ) barotropic velocities.
  44. ! vb )
  45. ! zeta surface elevation.
  46. ! zetado time derivative of *zeta*.
  47. ! psi barotropic stream function (pure diagnostic).
  48. !
  49. ! -------------------------------------------------------------------
  50. !
  51. ! Declaration of local constants/variables.
  52. ! -----------------------------------------
  53. !
  54. integer :: i,j,k
  55. !
  56. ! fields for predictor-corrector scheme
  57. !
  58. real (kind=8) :: ubold(ien,jen)
  59. real (kind=8) :: vbold(ien,jen)
  60. real (kind=8) :: zold(ien,jen)
  61. !
  62. !
  63. !* 1. Set initial constants and variables.
  64. ! ------------------------------------
  65. !
  66. zold(:,:) = zeta(:,:)
  67. !
  68. !* 2. Computation of normalized pressure.
  69. ! -----------------------------------
  70. !
  71. call press
  72. !
  73. !* 3. Computation of velocities.
  74. ! --------------------------
  75. !
  76. !* 3.2. *nsve*=2
  77. ! --------
  78. !
  79. ! Compute barotropic velocities.
  80. !
  81. call uvtrop2
  82. !
  83. ! Compute baroclinic velocities.
  84. !
  85. call geost2
  86. !
  87. ! predictor-corrector for sea ice/sea level
  88. !
  89. zeta(:,:) = zold(:,:)
  90. !
  91. call uvtrop2
  92. !
  93. ubold(:,:) = 0.0
  94. vbold(:,:) = 0.0
  95. do k=1,ken
  96. do j=3,jen-2
  97. do i=1,ien
  98. ubold(i,j)=ubold(i,j)+utot(i,j,k)*delta(i,j,k)
  99. vbold(i,j)=vbold(i,j)+vtot(i,j,k)*delta(i,j,k)
  100. end do
  101. end do
  102. end do
  103. do j=3,jen-2
  104. do i=1,ien
  105. if (depth(i,j)<0.5) cycle
  106. ubold(i,j)=ubold(i,j)/depth(i,j)
  107. vbold(i,j)=vbold(i,j)/depth(i,j)
  108. end do
  109. end do
  110. do k=1,ken
  111. do j=3,jen-2
  112. do i=1,ien
  113. utot(i,j,k)=(utot(i,j,k)+ub(i,j)-ubold(i,j))*wetvec(i,j,k)
  114. vtot(i,j,k)=(vtot(i,j,k)+vb(i,j)-vbold(i,j))*wetvec(i,j,k)
  115. end do
  116. end do
  117. end do
  118. !
  119. !
  120. !* 4. Compute vertical velocities.
  121. ! ----------------------------
  122. !
  123. call cont1
  124. return
  125. end subroutine diva