<p><b>dwj07@fsu.edu</b> 2011-10-25 14:37:08 -0600 (Tue, 25 Oct 2011)</p><p><br>
        -- BRANCH COMMIT --<br>
<br>
        Some small optimizations on the h_edge computation.<br>
<br>
        Removing the if statements around the 3rd and 4th order h_edge computations.<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/performance/src/core_ocean/mpas_ocn_tendency.F
===================================================================
--- branches/ocean_projects/performance/src/core_ocean/mpas_ocn_tendency.F        2011-10-25 19:25:17 UTC (rev 1134)
+++ branches/ocean_projects/performance/src/core_ocean/mpas_ocn_tendency.F        2011-10-25 20:37:08 UTC (rev 1135)
@@ -415,41 +415,31 @@
type (mesh_type), intent(in) :: grid
- integer :: iEdge, iCell, iVertex, k, cell1, cell2, vertex1, vertex2, eoe, i, j, cov
- real (kind=RKIND) :: flux, vorticity_abs, h_vertex, workpv, rho0Inv
+ integer :: iEdge, iCell, iVertex, k, cell1, cell2, vertex1, vertex2, eoe, i, j
+ integer :: boundaryMask, velMask, nCells, nEdges, nVertices, nVertLevels, vertexDegree, fCoef, err
- integer :: nCells, nEdges, nVertices, nVertLevels, vertexDegree, fCoef, err
+ integer, dimension(:), pointer :: nEdgesOnCell, nEdgesOnEdge, &
+ maxLevelCell, maxLevelEdgeTop, maxLevelEdgeBot, &
+ maxLevelVertexBot
+ integer, dimension(:,:), pointer :: cellsOnEdge, cellsOnVertex, &
+ verticesOnEdge, edgesOnEdge, edgesOnVertex,boundaryCell
+ real (kind=RKIND) :: d2fdx2_cell1, d2fdx2_cell2, coef_3rd_order, r_tmp, invAreaCell1, invAreaCell2, invAreaTri1, invAreaTri2, invLength, h_vertex
+ real (kind=RKIND), dimension(:), allocatable:: pTop
+
real (kind=RKIND), dimension(:), pointer :: &
- h_s, fVertex, fEdge, dvEdge, dcEdge, areaCell, areaTriangle, &
+ h_s, fVertex, dvEdge, dcEdge, areaCell, areaTriangle, &
hZLevel
real (kind=RKIND), dimension(:,:), pointer :: &
- weightsOnEdge, kiteAreasOnVertex, h_edge, h, u, v, w, pressure,&
- circulation, vorticity, ke, ke_edge, MontPot, wTop, &
- pv_edge, pv_vertex, pv_cell, gradPVn, gradPVt, divergence, &
- rho, temperature, salinity
- real (kind=RKIND), dimension(:,:,:), pointer :: tracers
- real (kind=RKIND), dimension(:), allocatable:: pTop
- real (kind=RKIND), dimension(:,:), allocatable:: div_u
- character :: c1*6
+ weightsOnEdge, kiteAreasOnVertex, h_edge, h, u, v, pressure,&
+ circulation, vorticity, ke, ke_edge, MontPot, &
+ pv_edge, pv_vertex, pv_cell, gradPVn, gradPVt, divergence, rho
+ real (kind=RKIND), dimension(:,:,:), pointer :: tracers, deriv_two
- integer, dimension(:,:), pointer :: cellsOnEdge, cellsOnVertex, &
- verticesOnEdge, edgesOnCell, edgesOnEdge, edgesOnVertex, &
- boundaryEdge, boundaryCell
- integer, dimension(:), pointer :: nEdgesOnCell, nEdgesOnEdge, &
- maxLevelCell, maxLevelEdgeTop, maxLevelEdgeBot, &
- maxLevelVertexBot, maxLevelVertexTop
- real (kind=RKIND) :: d2fdx2_cell1, d2fdx2_cell2
- real (kind=RKIND), dimension(:,:,:), pointer :: deriv_two
- real (kind=RKIND) :: coef_3rd_order
- real (kind=RKIND) :: r, h1, h2
- real (kind=RKIND) :: r_tmp, invAreaCell1, invAreaCell2, invAreaTri1, invAreaTri2, invLength
-
h => s % h % array
u => s % u % array
v => s % v % array
- wTop => s % wTop % array
h_edge => s % h_edge % array
circulation => s % circulation % array
vorticity => s % vorticity % array
@@ -462,7 +452,6 @@
gradPVn => s % gradPVn % array
gradPVt => s % gradPVt % array
rho => s % rho % array
- tracers => s % tracers % array
MontPot => s % MontPot % array
pressure => s % pressure % array
@@ -472,7 +461,6 @@
cellsOnVertex => grid % cellsOnVertex % array
verticesOnEdge => grid % verticesOnEdge % array
nEdgesOnCell => grid % nEdgesOnCell % array
- edgesOnCell => grid % edgesOnCell % array
nEdgesOnEdge => grid % nEdgesOnEdge % array
edgesOnEdge => grid % edgesOnEdge % array
edgesOnVertex => grid % edgesOnVertex % array
@@ -482,14 +470,12 @@
areaTriangle => grid % areaTriangle % array
h_s => grid % h_s % array
fVertex => grid % fVertex % array
- fEdge => grid % fEdge % array
hZLevel => grid % hZLevel % array
deriv_two => grid % deriv_two % array
maxLevelCell => grid % maxLevelCell % array
maxLevelEdgeTop => grid % maxLevelEdgeTop % array
maxLevelEdgeBot => grid % maxLevelEdgeBot % array
maxLevelVertexBot => grid % maxLevelVertexBot % array
- maxLevelVertexTop => grid % maxLevelVertexTop % array
nCells = grid % nCells
nEdges = grid % nEdges
@@ -497,7 +483,6 @@
nVertLevels = grid % nVertLevels
vertexDegree = grid % vertexDegree
- boundaryEdge => grid % boundaryEdge % array
boundaryCell => grid % boundaryCell % array
!
@@ -534,39 +519,27 @@
d2fdx2_cell1 = 0.0
d2fdx2_cell2 = 0.0
- !-- if not a boundary cell
- if(boundaryCell(k,cell1).eq.0.and.boundaryCell(k,cell2).eq.0) then
+ boundaryMask = .not.(boundaryCell(k,cell1) == 0 .and. boundaryCell(k,cell2) == 0)
- d2fdx2_cell1 = deriv_two(1,1,iEdge) * h(k,cell1)
- d2fdx2_cell2 = deriv_two(1,2,iEdge) * h(k,cell2)
+ d2fdx2_cell1 = deriv_two(1,1,iEdge) * h(k,cell1) * boundaryMask
+ d2fdx2_cell2 = deriv_two(1,2,iEdge) * h(k,cell2) * boundaryMask
- !-- all edges of cell 1
- do i=1, grid % nEdgesOnCell % array (cell1)
- d2fdx2_cell1 = d2fdx2_cell1 + &
- deriv_two(i+1,1,iEdge) * h(k,grid % CellsOnCell % array (i,cell1))
- end do
+ !-- all edges of cell 1
+ do i=1, nEdgesOnCell(cell1) * boundaryMask
+ d2fdx2_cell1 = d2fdx2_cell1 + &
+ deriv_two(i+1,1,iEdge) * h(k,grid % CellsOnCell % array (i,cell1))
+ end do
- !-- all edges of cell 2
- do i=1, grid % nEdgesOnCell % array (cell2)
- d2fdx2_cell2 = d2fdx2_cell2 + &
- deriv_two(i+1,2,iEdge) * h(k,grid % CellsOnCell % array (i,cell2))
- end do
+ !-- all edges of cell 2
+ do i=1, nEdgesOnCell(cell2) * boundaryMask
+ d2fdx2_cell2 = d2fdx2_cell2 + &
+ deriv_two(i+1,2,iEdge) * h(k,grid % CellsOnCell % array (i,cell2))
+ end do
- endif
+ velMask = 2*(u(k,iEdge) <= 0) - 1
- !-- if u > 0:
- if (u(k,iEdge) > 0) then
- h_edge(k,iEdge) = &
- 0.5*(h(k,cell1) + h(k,cell2)) &
- -(dcEdge(iEdge) **2) * (d2fdx2_cell1 + d2fdx2_cell2) / 12. &
- -(dcEdge(iEdge) **2) * coef_3rd_order*(d2fdx2_cell1 - d2fdx2_cell2) / 12.
- !-- else u <= 0:
- else
- h_edge(k,iEdge) = &
- 0.5*(h(k,cell1) + h(k,cell2)) &
- -(dcEdge(iEdge) **2) * (d2fdx2_cell1 + d2fdx2_cell2) / 12. &
- +(dcEdge(iEdge) **2) * coef_3rd_order*(d2fdx2_cell1 - d2fdx2_cell2) / 12.
- end if
+ h_edge(k,iEdge) = 0.5*(h(k,cell1) + h(k,cell2)) - (dcEdge(iEdge) **2) * (d2fdx2_cell1 + d2fdx2_cell2) / 12. &
+ + velMask * (dcEdge(iEdge) **2) * coef_3rd_order*(d2fdx2_cell1 - d2fdx2_cell2) / 12.
end do ! do k
end do ! do iEdge
@@ -582,26 +555,23 @@
d2fdx2_cell1 = 0.0
d2fdx2_cell2 = 0.0
- !-- if not a boundary cell
- if(boundaryCell(k,cell1).eq.0.and.boundaryCell(k,cell2).eq.0) then
+ boundaryMask = .not.(boundaryCell(k,cell1) == 0 .and. boundaryCell(k,cell2) == 0)
- d2fdx2_cell1 = deriv_two(1,1,iEdge) * h(k,cell1)
- d2fdx2_cell2 = deriv_two(1,2,iEdge) * h(k,cell2)
+ d2fdx2_cell1 = deriv_two(1,1,iEdge) * h(k,cell1) * boundaryMask
+ d2fdx2_cell2 = deriv_two(1,2,iEdge) * h(k,cell2) * boundaryMask
- !-- all edges of cell 1
- do i=1, grid % nEdgesOnCell % array (cell1)
- d2fdx2_cell1 = d2fdx2_cell1 + &
- deriv_two(i+1,1,iEdge) * h(k,grid % CellsOnCell % array (i,cell1))
- end do
+ !-- all edges of cell 1
+ do i=1, nEdgesOnCell(cell1) * boundaryMask
+ d2fdx2_cell1 = d2fdx2_cell1 + &
+ deriv_two(i+1,1,iEdge) * h(k,grid % CellsOnCell % array (i,cell1))
+ end do
- !-- all edges of cell 2
- do i=1, grid % nEdgesOnCell % array (cell2)
- d2fdx2_cell2 = d2fdx2_cell2 + &
- deriv_two(i+1,2,iEdge) * h(k,grid % CellsOnCell % array (i,cell2))
- end do
+ !-- all edges of cell 2
+ do i=1, nEdgesOnCell(cell2) * boundaryMask
+ d2fdx2_cell2 = d2fdx2_cell2 + &
+ deriv_two(i+1,2,iEdge) * h(k,grid % CellsOnCell % array (i,cell2))
+ end do
- endif
-
h_edge(k,iEdge) = &
0.5*(h(k,cell1) + h(k,cell2)) &
-(dcEdge(iEdge) **2) * (d2fdx2_cell1 + d2fdx2_cell2) / 12.
</font>
</pre>