Skip to content

Commit

Permalink
Use alternate weighted z-star formulation
Browse files Browse the repository at this point in the history
  • Loading branch information
cbegeman committed May 23, 2024
1 parent aaa852a commit f9d2573
Showing 1 changed file with 25 additions and 14 deletions.
39 changes: 25 additions & 14 deletions components/mpas-ocean/src/shared/mpas_ocn_thick_ale.F
Original file line number Diff line number Diff line change
Expand Up @@ -125,10 +125,12 @@ subroutine ocn_ALE_thickness(verticalMeshPool, SSH, ALE_thickness, &
nCells ! number of cells

real (kind=RKIND) :: &
weightSum, &! sum of weights in vertical
thicknessSum, &! total thickness
remainder, &! track remainder in mix/max alteration
newThickness ! temp used during min/max adjustment
weightSum, &! sum of weights in vertical over 1 to maxLevelCell
activeWeightSum, &! sum of weights in vertical over minLevelCell to maxLevelCell
restingThicknessSum,&! total resting thickness
newThicknessSum, &! total thickness at current time
remainder, &! track remainder in mix/max alteration
newThickness ! temp used during min/max adjustment

real (kind=RKIND), dimension(:), allocatable :: &
prelim_ALE_thickness, & ! ALE thickness at new time
Expand Down Expand Up @@ -169,34 +171,43 @@ subroutine ocn_ALE_thickness(verticalMeshPool, SSH, ALE_thickness, &

#ifdef MPAS_OPENACC
!xacc parallel loop &
!xacc present(ALE_thickness, SSH, restingThickness, &
!xacc present(ALE_thickness, SSH, restingThicknessSum, &
!xacc newThicknessSum, weightSum, activeWeightSum, &
!xacc minLevelCell, maxLevelCell, &
!xacc vertCoordMovementWeights) &
!xacc private(k, kMin, kMax, thicknessSum)
!xacc private(k, kMin, kMax, restingThicknessSum, &
!xacc newThicknessSum, weightSum, activeWeightSum)
#else
!$omp parallel
!$omp do schedule(runtime) &
!$omp private(k, kMin, kMax, thicknessSum)
!$omp private(k, kMin, kMax, restingThicknessSum, &
!$omp newThicknessSum, weightSum, activeWeightSum)
#endif
do iCell = 1, nCells
kMax = maxLevelCell(iCell)
kMin = minLevelCell(iCell)

thicknessSum = 1e-14_RKIND
restingThicknessSum = 1e-14_RKIND
weightSum = 0.0_RKIND
activeWeightSum = 1e-14_RKIND
do k = kMin, kMax
thicknessSum = thicknessSum &
+ vertCoordMovementWeights(k) &
* restingThickness(k,iCell)
restingThicknessSum = restingThicknessSum &
+ restingThickness(k,iCell)
activeWeightSum = activeWeightSum + vertCoordMovementWeights(k)
end do
do k = 1, kMax
weightSum = weightSum + vertCoordMovementWeights(k)
end do
newThicknessSum = bottomDepth(iCell) + ssh(iCell)

! Note that restingThickness is nonzero, and remaining
! terms are perturbations about zero.
! This is equation 4 and 6 in Petersen et al 2015,
! but with eqn 6
do k = kMin, kMax
ALE_thickness(k,iCell) = restingThickness(k,iCell) &
+ (SSH(iCell)*vertCoordMovementWeights(k)* &
restingThickness(k,iCell) )/thicknessSum
ALE_thickness(k, iCell) = (newThicknessSum / restingThicknessSum ) &
* (weightSum / activeWeightSum) &
* vertCoordMovementWeights(k) * restingThickness(k, iCell)
end do
enddo
#ifndef MPAS_OPENACC
Expand Down

0 comments on commit f9d2573

Please sign in to comment.