<p><b>duda</b> 2010-04-28 17:11:14 -0600 (Wed, 28 Apr 2010)</p><p>BRANCH COMMIT<br>
<br>
Fix bug in cam_restart_to_mpas() interface routine: we had been calling<br>
cam_inidat_to_mpas() from cam_restart_to_mpas() to initialize all but the<br>
omega field, but passed potential temperature rather than temperature to<br>
cam_inidat_to_mpas(). Art has agreed to modify the interface specification<br>
so that both cam_restart_to_mpas() and cam_inidat_to_mpas() are given<br>
potential temperature, so we just need to update cam_inidat_to_mpas()<br>
to use potential temperature.<br>
<br>
Also, modify the code in the driver that writes out test input files so that<br>
these files contain potential temperature rather than temperature.<br>
<br>
M driver/mpas.F<br>
M driver_cam_interface/module_mpas_cam_interface.F<br>
</p><hr noshade><pre><font color="gray">Modified: branches/mpas_cam_coupling/src/driver/mpas.F
===================================================================
--- branches/mpas_cam_coupling/src/driver/mpas.F        2010-04-28 21:45:59 UTC (rev 220)
+++ branches/mpas_cam_coupling/src/driver/mpas.F        2010-04-28 23:11:14 UTC (rev 221)
@@ -60,10 +60,10 @@
do k=1, domain % blocklist % mesh % nVertLevels
T(iCell,k) = domain % blocklist % time_levs(1) % state % theta % array(k,iCell)
- T(iCell,k) = T(iCell,k) * (0.5 * (domain % blocklist % time_levs(1) % state % pressure % array(k, iCell) + &
- domain % blocklist % time_levs(1) % state % pressure % array(k+1, iCell)) &
- / 100000.0 &
- ) ** (287.0 / 1003.0)
+! T(iCell,k) = T(iCell,k) * (0.5 * (domain % blocklist % time_levs(1) % state % pressure % array(k, iCell) + &
+! domain % blocklist % time_levs(1) % state % pressure % array(k+1, iCell)) &
+! / 100000.0 &
+! ) ** (287.0 / 1003.0)
do j=1,domain % blocklist % mesh % nEdgesOnCell % array(iCell)
iEdge = domain % blocklist % mesh % edgesOnCell % array(j,iCell)
!NOTE: Might need to adjust sign +/-
Modified: branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F
===================================================================
--- branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F        2010-04-28 21:45:59 UTC (rev 220)
+++ branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F        2010-04-28 23:11:14 UTC (rev 221)
@@ -127,7 +127,7 @@
subroutine cam_inidat_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Pref, Ptop, &
- Ae, Be, Ac, Bc, Psd, Phis, T, U, Tracer)
+ Ae, Be, Ac, Bc, Psd, Phis, Theta, U, Tracer)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! SUBROUTINE MPAS_INIDAT_TO_MPAS
!
@@ -141,7 +141,7 @@
! Ac, Bc - A and B arrays at vertical centers
! Psd - dry surface pressure (Pa)
! Phis - surface geopotential (m^2/s^2)
- ! T - temperature (K)
+ ! Theta - potential temperature (K)
! U - normal velocity at edges (m/s)
! Tracer - tracer mixing ratios (kg/kg)
!
@@ -159,7 +159,7 @@
real (kind=RKIND), dimension(Plev), intent(in) :: Ac, Bc
real (kind=RKIND), dimension(Numcols), intent(in) :: Psd
real (kind=RKIND), dimension(Numcols), intent(in) :: Phis
- real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: T
+ real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: Theta
real (kind=RKIND), dimension(Numcols,MaxEdges,Plev), intent(in) :: U
real (kind=RKIND), dimension(Numcols,Plev,Pcnst), intent(in) :: Tracer
@@ -169,7 +169,7 @@
real (kind=RKIND), dimension(Plev+1) :: znuc, znuv, bn, divh, dpn
type (block_type), pointer :: block
real (kind=RKIND), dimension(:), pointer :: surface_pressure, rdnu, rdnw, fnm, fnp, dbn, dnu, dnw
- real (kind=RKIND), dimension(:,:), pointer :: theta, u_normal, pressure, geopotential, alpha, h
+ real (kind=RKIND), dimension(:,:), pointer :: t_pot, u_normal, pressure, geopotential, alpha, h
real (kind=RKIND), dimension(:,:,:), pointer :: scalars
write(0,*) 'Called CAM_INIDAT_TO_MPAS'
@@ -177,7 +177,7 @@
block => domain % blocklist
surface_pressure => block % time_levs(1) % state % surface_pressure % array
- theta => block % time_levs(1) % state % theta % array
+ t_pot => block % time_levs(1) % state % theta % array
u_normal => block % time_levs(1) % state % u % array
pressure => block % time_levs(1) % state % pressure % array
geopotential => block % time_levs(1) % state % geopotential % array
@@ -221,7 +221,7 @@
write(0,*) 'min/max for Bc = ', minval(Bc), maxval(Bc)
write(0,*) 'min/max for Psd = ', minval(Psd), maxval(Psd)
write(0,*) 'min/max for Phis = ', minval(Phis), maxval(Phis)
- write(0,*) 'min/max for T = ', minval(T), maxval(T)
+ write(0,*) 'min/max for Theta = ', minval(Theta), maxval(Theta)
write(0,*) 'min/max for U = ', minval(U), maxval(U)
write(0,*) 'min/max for Tracer = ', minval(Tracer), maxval(Tracer)
@@ -308,7 +308,7 @@
!
h(:,:) = 0.0
pressure(:,:) = 0.0
- theta(:,:) = 0.0
+ t_pot(:,:) = 0.0
alpha(:,:) = 0.0
geopotential(:,:) = 0.0
do iCell=1,Numcols
@@ -322,8 +322,8 @@
end do
do k=1,Plev
- theta(k,iCell) = T(iCell,k) * (0.5*(pressure(k,iCell) + pressure(k+1,iCell))/p0)**(-rgas/cp)
- alpha(k,iCell) = (rgas/p0) * theta(k,iCell) * (0.5*(pressure(k,iCell)+pressure(k+1,iCell))/p0)**cvpm
+ t_pot(k,iCell) = Theta(iCell,k)
+ alpha(k,iCell) = (rgas/p0) * Theta(iCell,k) * (0.5*(pressure(k,iCell)+pressure(k+1,iCell))/p0)**cvpm
end do
geopotential(1,iCell) = Phis(iCell)
@@ -376,9 +376,9 @@
subroutine cam_restart_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Pref, Ptop, &
- Ae, Be, Ac, Bc, Psd, Phis, T, U, W, Tracer)
+ Ae, Be, Ac, Bc, Psd, Phis, Theta, U, W, Tracer)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
- ! SUBROUTINE MPAS_RESTART_TO_MPAS
+ ! SUBROUTINE CAM_RESTART_TO_MPAS
!
! Input: Numcols - number of columns/cells
! MaxEdges - maximum number of edges per cell
@@ -390,7 +390,7 @@
! Ac, Bc - A and B arrays at vertical centers
! Psd - dry surface pressure (Pa)
! Phis - surface geopotential (m^2/s^2)
- ! T - temperature (K)
+ ! Theta - potential temperature (K)
! U - normal velocity at edges (m/s)
! W - vertical velocity (Pa/s)
! Tracer - tracer mixing ratios (kg/kg)
@@ -408,7 +408,7 @@
real (kind=RKIND), dimension(Plev), intent(in) :: Ac, Bc
real (kind=RKIND), dimension(Numcols), intent(in) :: Psd
real (kind=RKIND), dimension(Numcols), intent(in) :: Phis
- real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: T
+ real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: Theta
real (kind=RKIND), dimension(Numcols,MaxEdges,Plev), intent(in) :: U
real (kind=RKIND), dimension(Numcols,Plev+1), intent(in) :: W
real (kind=RKIND), dimension(Numcols,Plev,Pcnst), intent(in) :: Tracer
@@ -442,7 +442,7 @@
! Initialize most variables using cam_inidat_to_mpas
!
call cam_inidat_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Pref, Ptop, &
- Ae, Be, Ac, Bc, Psd, Phis, T, U, Tracer)
+ Ae, Be, Ac, Bc, Psd, Phis, Theta, U, Tracer)
do iCell=1,Numcols
do k=1,Plev+1
</font>
</pre>