<p><b>mpetersen@lanl.gov</b> 2010-04-22 08:51:06 -0600 (Thu, 22 Apr 2010)</p><p>Add T,S, test tracer variables, and linear equation of state<br>
</p><hr noshade><pre><font color="gray">Modified: branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/Registry
===================================================================
--- branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/Registry        2010-04-21 23:30:52 UTC (rev 200)
+++ branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/Registry        2010-04-22 14:51:06 UTC (rev 201)
@@ -78,7 +78,6 @@
var real fEdge ( nEdges ) iro fEdge - -
var real fVertex ( nVertices ) iro fVertex - -
var real h_s ( nCells ) iro h_s - -
-var real rho ( nVertLevels nCells Time ) iro rho - -
# Arrays required for reconstruction of velocity field
var real coeffs_reconstruct ( R3 maxEdges nCells ) - coeffs_reconstruct - -
@@ -90,6 +89,7 @@
# Prognostic variables: read from input, saved in restart, and written to output
var real u ( nVertLevels nEdges Time ) iro u - -
var real h ( nVertLevels nCells Time ) iro h - -
+var real rho ( nVertLevels nCells Time ) io rho - -
var real temperature ( nVertLevels nCells Time ) iro temperature tracers dynamics
var real salinity ( nVertLevels nCells Time ) iro salinity tracers dynamics
var real tracer1 ( nVertLevels nCells Time ) iro tracer1 tracers testing
Modified: branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_global_diagnostics.F
===================================================================
--- branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_global_diagnostics.F        2010-04-21 23:30:52 UTC (rev 200)
+++ branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_global_diagnostics.F        2010-04-22 14:51:06 UTC (rev 201)
@@ -36,7 +36,7 @@
real (kind=RKIND) :: areaCellGlobal, areaEdgeGlobal, areaTriangleGlobal
real (kind=RKIND), dimension(:), pointer :: areaCell, dcEdge, dvEdge, areaTriangle, areaEdge
real (kind=RKIND), dimension(:,:), pointer :: h, u, v, h_edge, circulation, vorticity, ke, pv_edge, pv_vertex, &
- pv_cell, gradPVn, gradPVt, zmid, zbot, pmid, pbot, MontPot, w
+ pv_cell, gradPVn, gradPVt, zmid, zbot, pmid, pbot, MontPot, w, rho
real (kind=RKIND), dimension(:,:,:), pointer :: tracers
real (kind=RKIND) :: volumeCellGlobal, volumeEdgeGlobal, CFLNumberGlobal
real (kind=RKIND) :: localCFL, localSum
@@ -63,6 +63,7 @@
h => state % h % array
u => state % u % array
+ rho => state % rho % array
tracers => state % tracers % array
v => state % v % array
w => state % w % array
@@ -361,6 +362,8 @@
if (nVertLevels.gt.1) &
print '(a,10es25.15)', 'k=2+ min h max h',&
minval(h(2:nVertLevels,1:nCellsSolve)),maxval(h(2:nVertLevels,1:nCellsSolve))
+ print '(a,10es25.15)', ' min rho max rho',&
+ minval(rho(:,1:nCellsSolve)),maxval(rho(:,1:nCellsSolve))
do k=1,4
print '(a,i5,10es25.15)', 'tracer min max',&
Modified: branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_test_cases.F
===================================================================
--- branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_test_cases.F        2010-04-21 23:30:52 UTC (rev 200)
+++ branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_test_cases.F        2010-04-22 14:51:06 UTC (rev 201)
@@ -78,7 +78,7 @@
stop
end if
-!print *, temperature_INDEX
+print *, 'index_temperature',index_temperature
! mrp 100121:
!
@@ -116,10 +116,16 @@
tracers = 0.0
do iCell = 1,nCells
do iLevel = 1,nVertLevels
- tracers(3,iLevel,iCell) = xCell(iCell)
- if (xCell(iCell)<500*1e3) tracers(4,iLevel,iCell) = 1.0
- tracers(1,iLevel,iCell) = xCell(iCell)
- tracers(2,iLevel,iCell) = 1.0
+ !tracers(1,iLevel,iCell) = 2.0
+ !tracers(2,iLevel,iCell) = 1.0
+ !tracers(3,iLevel,iCell) = xCell(iCell)
+ !if (xCell(iCell)<500*1e3) tracers(4,iLevel,iCell) = 1.0
+! tracers(2,iLevel,iCell) = 5.0 ! temperature
+! tracers(3,iLevel,iCell) = 35.0 ! salinity
+ tracers(index_temperature,iLevel,iCell) = &
+ 10.0* xCell(iCell)/2500e3 ! temperature
+ tracers(index_salinity,iLevel,iCell) = &
+ 34.0 + 2* yCell(iCell)/4000e3! salinity
enddo
enddo
Modified: branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_time_integration.F
===================================================================
--- branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_time_integration.F        2010-04-21 23:30:52 UTC (rev 200)
+++ branches/ocean_projects/z_level_mrp/mpas/src/core_ocean/module_time_integration.F        2010-04-22 14:51:06 UTC (rev 201)
@@ -543,18 +543,20 @@
tend_tr(:,:,:) = 0.0
do iEdge=1,nEdges
- cell1 = cellsOnEdge(1,iEdge)
- cell2 = cellsOnEdge(2,iEdge)
- if (cell1 > 0 .and. cell2 > 0) then
- do k=1,nVertLevels
- do iTracer=1,num_tracers
- tracer_edge = 0.5 * (tracers(iTracer,k,cell1) + tracers(iTracer,k,cell2))
- flux = u(k,iEdge) * dvEdge(iEdge) * h_edge(k,iEdge) * tracer_edge
- tend_tr(iTracer,k,cell1) = tend_tr(iTracer,k,cell1) - flux
- tend_tr(iTracer,k,cell2) = tend_tr(iTracer,k,cell2) + flux
- end do
+ cell1 = cellsOnEdge(1,iEdge)
+ cell2 = cellsOnEdge(2,iEdge)
+ if (cell1 > 0 .and. cell2 > 0) then
+ do k=1,nVertLevels
+ do iTracer=1,num_tracers
+ tracer_edge = 0.5 * ( tracers(iTracer,k,cell1) &
+ + tracers(iTracer,k,cell2))
+ flux = u(k,iEdge) * dvEdge(iEdge) * h_edge(k,iEdge) &
+ * tracer_edge
+ tend_tr(iTracer,k,cell1) = tend_tr(iTracer,k,cell1) - flux
+ tend_tr(iTracer,k,cell2) = tend_tr(iTracer,k,cell2) + flux
end do
- end if
+ end do
+ end if
end do
! mrp 100409 computation of h tendency was, only had div.(huT) term
@@ -641,7 +643,7 @@
vh, weightsOnEdge, kiteAreasOnVertex, h_edge, h, u, v, w, &
tend_h, tend_u, circulation, vorticity, ke, &
pv_edge, pv_vertex, pv_cell, gradPVn, gradPVt, divergence, &
- zmid, zbot, pmid, pbot, MontPot, rho
+ zmid, zbot, pmid, pbot, MontPot, rho, temperature, salinity
character :: c1*6
integer, dimension(:,:), pointer :: cellsOnEdge, cellsOnVertex, verticesOnEdge, edgesOnCell, edgesOnEdge, edgesOnVertex, uBC
@@ -668,6 +670,8 @@
gradPVt => s % gradPVt % array
!mrp 100112:
rho => s % rho % array
+ temperature => s % tracers % array(index_temperature,:,:)
+ salinity => s % tracers % array(index_salinity,:,:)
zmid => s % zmid % array
zbot => s % zbot % array
pmid => s % pmid % array
@@ -943,6 +947,21 @@
enddo
endif
+ !mrp 100324:
+ !
+ ! equation of state
+ !
+ do iCell=1,nCells
+ do k=1,nVertLevels
+ ! Linear equation of state, for the time being
+ rho(k,iCell) = 1000.0*( 1.0 &
+ - 2.5e-4*temperature(k,iCell) &
+ + 7.6e-4*salinity(k,iCell))
+ end do
+ end do
+ !mrp 100324 end
+
+
!mrp 100112:
!
! Compute mid- and bottom-depth of each layer, from bottom up.
</font>
</pre>