<p><b>ringler@lanl.gov</b> 2011-05-18 11:10:25 -0600 (Wed, 18 May 2011)</p><p><br>
add option to scale del2 and del4 viscosity based on mesh grid spacing<br>
<br>
the option is turned on with "config_h_ScaleWithMesh = .true." in namelist (default is false)<br>
<br>
the scaling is done based on meshDensity that is now a part of Registry. This field needs to be in the grid.nc file for the option to work.<br>
<br>
the current scaling changes the del4 parameter by a factor of 10 for each doubling of grid spacing.<br>
<br>
note, when using this option the :config_h_mom_eddy_visc2" and "config_h_mom_eddy_visc4" parameters should be set to appropriately for the finest region of the mesh. The scaling will increase these parameters in the coarser regions.<br>
</p><hr noshade><pre><font color="gray">Modified: trunk/mpas/namelist.input.sw
===================================================================
--- trunk/mpas/namelist.input.sw        2011-05-18 15:02:58 UTC (rev 840)
+++ trunk/mpas/namelist.input.sw        2011-05-18 17:10:25 UTC (rev 841)
@@ -5,6 +5,7 @@
config_ntimesteps = 7500
config_output_interval = 500
config_stats_interval = 0
+ config_h_ScaleWithMesh = .false.
config_h_mom_eddy_visc2 = 0.0
config_h_mom_eddy_visc4 = 0.0
config_h_tracer_eddy_diff2 = 0.0
Modified: trunk/mpas/src/core_sw/Registry
===================================================================
--- trunk/mpas/src/core_sw/Registry        2011-05-18 15:02:58 UTC (rev 840)
+++ trunk/mpas/src/core_sw/Registry        2011-05-18 17:10:25 UTC (rev 841)
@@ -7,6 +7,7 @@
namelist integer sw_model config_ntimesteps 7500
namelist integer sw_model config_output_interval 500
namelist integer sw_model config_stats_interval 100
+namelist logical sw_model config_h_ScaleWithMesh false
namelist real sw_model config_h_mom_eddy_visc2 0.0
namelist real sw_model config_h_mom_eddy_visc4 0.0
namelist real sw_model config_h_tracer_eddy_diff2 0.0
@@ -52,6 +53,7 @@
var persistent real xCell ( nCells ) 0 iro xCell mesh - -
var persistent real yCell ( nCells ) 0 iro yCell mesh - -
var persistent real zCell ( nCells ) 0 iro zCell mesh - -
+var persistent real meshDensity ( nCells ) 0 iro meshDensity mesh - -
var persistent integer indexToCellID ( nCells ) 0 iro indexToCellID mesh - -
var persistent real latEdge ( nEdges ) 0 iro latEdge mesh - -
Modified: trunk/mpas/src/core_sw/module_time_integration.F
===================================================================
--- trunk/mpas/src/core_sw/module_time_integration.F        2011-05-18 15:02:58 UTC (rev 840)
+++ trunk/mpas/src/core_sw/module_time_integration.F        2011-05-18 17:10:25 UTC (rev 841)
@@ -258,7 +258,7 @@
real (kind=RKIND) :: flux, vorticity_abs, workpv, q, upstream_bias
integer :: nCells, nEdges, nVertices, nVertLevels
- real (kind=RKIND), dimension(:), pointer :: h_s, fVertex, fEdge, dvEdge, dcEdge, areaCell, areaTriangle
+ real (kind=RKIND), dimension(:), pointer :: h_s, fVertex, fEdge, dvEdge, dcEdge, areaCell, areaTriangle, meshDensity
real (kind=RKIND), dimension(:,:), pointer :: vh, weightsOnEdge, kiteAreasOnVertex, h_edge, h, u, v, tend_h, tend_u, &
circulation, vorticity, ke, pv_edge, divergence, h_vertex
integer, dimension(:,:), pointer :: cellsOnEdge, cellsOnVertex, verticesOnEdge, edgesOnCell, edgesOnEdge, edgesOnVertex
@@ -273,6 +273,8 @@
real (kind=RKIND), parameter :: rho_ref = 1000.0
real (kind=RKIND) :: ke_edge
+ real (kind=RKIND), allocatable, dimension(:) :: meshScalingDel2, meshScalingDel4
+
h => s % h % array
u => s % u % array
v => s % v % array
@@ -312,7 +314,25 @@
u_src => grid % u_src % array
+ meshDensity => grid % meshDensity % array
+
!
+ ! Compute the scaling factors to be used in the del2 and del4 dissipation
+ !
+ allocate(meshScalingDel2(nEdges))
+ allocate(meshScalingDel4(nEdges))
+ meshScalingDel2(:) = 1.0
+ meshScalingDel4(:) = 1.0
+ if(config_h_ScaleWithMesh) then
+ do iEdge=1,nEdges
+ cell1 = cellsOnEdge(1,iEdge)
+ cell2 = cellsOnEdge(2,iEdge)
+ meshScalingDel2(iEdge) = ( (meshDensity(cell1) + meshDensity(cell2) )/2.0)**(5.0/12.0)
+ meshScalingDel4(iEdge) = (meshScalingDel2(iEdge))**2
+ enddo
+ endif
+
+ !
! Compute height tendency for each cell
!
tend_h(:,:) = 0.0
@@ -399,7 +419,7 @@
do k=1,nVertLevels
u_diffusion = ( divergence(k,cell2) - divergence(k,cell1) ) / dcEdge(iEdge) &
-(vorticity(k,vertex2) - vorticity(k,vertex1) ) / dvEdge(iEdge)
- u_diffusion = config_h_mom_eddy_visc2 * u_diffusion
+ u_diffusion = meshScalingDel2(iEdge) * config_h_mom_eddy_visc2 * u_diffusion
tend_u(k,iEdge) = tend_u(k,iEdge) + u_diffusion
end do
end do
@@ -487,7 +507,7 @@
-( delsq_vorticity(k,vertex2) &
- delsq_vorticity(k,vertex1) ) / dvEdge(iEdge)
- u_diffusion = config_h_mom_eddy_visc4 * u_diffusion
+ u_diffusion = meshScalingDel4(iEdge) * config_h_mom_eddy_visc4 * u_diffusion
tend_u(k,iEdge) = tend_u(k,iEdge) - u_diffusion
end do
@@ -521,6 +541,9 @@
*sqrt(2.0*ke_edge)/h_edge(1,iEdge)
end do
endif
+
+ deallocate(meshScalingDel2)
+ deallocate(meshScalingDel4)
end subroutine compute_tend
</font>
</pre>