<p><b>duda</b> 2010-03-18 11:47:52 -0600 (Thu, 18 Mar 2010)</p><p>BRANCH COMMIT<br>
<br>
Fill in details of MPAS-CAM interface routine arguments based on <br>
LLNL's interface document. Update test driver program accordingly.<br>
<br>
M src/driver_cam_interface/module_mpas_cam_interface.F<br>
M src/driver_cam_interface/test_driver.F<br>
</p><hr noshade><pre><font color="gray">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-03-17 20:33:00 UTC (rev 144)
+++ branches/mpas_cam_coupling/src/driver_cam_interface/module_mpas_cam_interface.F        2010-03-18 17:47:52 UTC (rev 145)
@@ -12,6 +12,9 @@
real (kind=RKIND) :: dt
integer :: itimestep
+ ! TEMPORARY: Once we are initializing MPAS fields from CAM-supplied
+ ! arrays, we won't need to call mpas_setup_test_case to get initial
+ ! fields
interface mpas_setup_test_case
subroutine mpas_setup_test_case(domain)
use grid_types
@@ -49,14 +52,16 @@
call allocate_domain(domain, dminfo)
-! NOTE: Really, we only want to read in the grid information, since any state
-! variables will be overwritten by cam_inidat_to_mpas() or cam_restart_to_mpas()
+ ! NOTE: Really, we only want to read in the grid information, since any state
+ ! variables will be overwritten by cam_inidat_to_mpas() or cam_restart_to_mpas();
+ ! however, with a suitable registry, this shouldn't be a problem
call input_state_for_domain(domain)
end subroutine mpas_init1
- subroutine cam_inidat_to_mpas()
+ subroutine cam_inidat_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Pref, Ptop, &
+ Ae, Be, Ac, Bc, Psd, Phis, T, U, Tracer)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! SUBROUTINE MPAS_INIDAT_TO_MPAS
!
@@ -77,17 +82,76 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
implicit none
+
+ integer, intent(in) :: Numcols
+ integer, intent(in) :: MaxEdges
+ integer, intent(in) :: Plev
+ integer, intent(in) :: Pcnst
+ real (kind=RKIND), intent(in) :: Pref
+ real (kind=RKIND), intent(in) :: Ptop
+ real (kind=RKIND), dimension(Plev+1), intent(in) :: Ae, Be
+ 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,MaxEdges,Plev), intent(in) :: U
+ real (kind=RKIND), dimension(Numcols,Plev,Pcnst), intent(in) :: Tracer
+
+ type (block_type), pointer :: block
write(0,*) 'Called CAM_INIDAT_TO_MPAS'
+ block => domain % blocklist
+
+ !
+ ! Perform basic sanity check on expected and available field dimensions
+ !
+ if (Numcols /= block % mesh % nCellsSolve) then
+ write(0,*) 'Error: mismatch between Numcols and nCellsSolve: ', Numcols, block % mesh % nCellsSolve
+ end if
+ if (Plev /= block % mesh % nVertLevels) then
+ write(0,*) 'Error: mismatch between Plev and nVertLevels: ', Plev, block % mesh % nVertLevels
+ end if
+ if (Pcnst /= num_scalars) then
+ write(0,*) 'Error: mismatch between Pcnst and num_scalars: ', Pcnst, num_scalars
+ end if
+
+
+ !
+ ! Fill in MPAS fields in block % time_levs(1) % state arrays from arguments
+ !
+
! TEMPORARY: until we are actually given fields, just call one of the test case
! setup routines to initialize the MPAS prognostic fields
call mpas_setup_test_case(domain)
+
+
+ !
+ ! Do a ghost-cell update
+ !
+ call dmpar_exch_halo_field1dReal(domain % dminfo, block % time_levs(1) % state % surface_pressure % array(:), &
+ block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % theta % array(:,:), &
+ block % mesh % nVertLevels, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % u % array(:,:), &
+ block % mesh % nVertLevels, block % mesh % nEdges, &
+ block % parinfo % edgesToSend, block % parinfo % edgesToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % pressure % array(:,:), &
+ block % mesh % nVertLevels+1, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % geopotential % array(:,:), &
+ block % mesh % nVertLevels+1, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field3dReal(domain % dminfo, block % time_levs(1) % state % scalars % array(:,:,:), &
+ num_scalars, block % mesh % nVertLevels, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
end subroutine cam_inidat_to_mpas
- subroutine cam_restart_to_mpas()
+ subroutine cam_restart_to_mpas(Numcols, MaxEdges, Plev, Pcnst, Psd, Phis, T, U, W, Tracer)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! SUBROUTINE MPAS_RESTART_TO_MPAS
!
@@ -104,12 +168,71 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
implicit none
+
+ integer, intent(in) :: Numcols
+ integer, intent(in) :: MaxEdges
+ integer, intent(in) :: Plev
+ integer, intent(in) :: Pcnst
+ 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,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
+
+ type (block_type), pointer :: block
write(0,*) 'Called CAM_RESTART_TO_MPAS'
+ block => domain % blocklist
+
+ !
+ ! Perform basic sanity check on expected and available field dimensions
+ !
+ if (Numcols /= block % mesh % nCellsSolve) then
+ write(0,*) 'Error: mismatch between Numcols and nCellsSolve: ', Numcols, block % mesh % nCellsSolve
+ end if
+ if (Plev /= block % mesh % nVertLevels) then
+ write(0,*) 'Error: mismatch between Plev and nVertLevels: ', Plev, block % mesh % nVertLevels
+ end if
+ if (Pcnst /= num_scalars) then
+ write(0,*) 'Error: mismatch between Pcnst and num_scalars: ', Pcnst, num_scalars
+ end if
+
+
+ !
+ ! Fill in MPAS fields in block % time_levs(1) % state arrays from arguments
+ !
+
! TEMPORARY: until we are actually given fields, just call one of the test case
! setup routines to initialize the MPAS prognostic fields
call mpas_setup_test_case(domain)
+
+
+ !
+ ! Do a ghost-cell update
+ !
+ call dmpar_exch_halo_field1dReal(domain % dminfo, block % time_levs(1) % state % surface_pressure % array(:), &
+ block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % theta % array(:,:), &
+ block % mesh % nVertLevels, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % u % array(:,:), &
+ block % mesh % nVertLevels, block % mesh % nEdges, &
+ block % parinfo % edgesToSend, block % parinfo % edgesToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % pressure % array(:,:), &
+ block % mesh % nVertLevels+1, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % geopotential % array(:,:), &
+ block % mesh % nVertLevels+1, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field2dReal(domain % dminfo, block % time_levs(1) % state % ww % array(:,:), &
+ block % mesh % nVertLevels+1, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
+ call dmpar_exch_halo_field3dReal(domain % dminfo, block % time_levs(1) % state % scalars % array(:,:,:), &
+ num_scalars, block % mesh % nVertLevels, block % mesh % nCells, &
+ block % parinfo % cellsToSend, block % parinfo % cellsToRecv)
end subroutine cam_restart_to_mpas
@@ -141,7 +264,7 @@
end subroutine mpas_init2
- subroutine mpas_to_cam()
+ subroutine mpas_to_cam(Numcols, Plev, Pcnst, Psd, Phis, T, Ux, Uy, Omega, Tracer)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! SUBROUTINE MPAS_TO_CAM
!
@@ -161,12 +284,52 @@
implicit none
+ integer, intent(in) :: Numcols
+ integer, intent(in) :: Plev
+ integer, intent(in) :: Pcnst
+ real (kind=RKIND), dimension(Numcols), intent(out) :: Psd
+ real (kind=RKIND), dimension(Numcols), intent(out) :: Phis
+ real (kind=RKIND), dimension(Numcols,Plev), intent(out) :: T
+ real (kind=RKIND), dimension(Numcols,Plev), intent(out) :: Ux
+ real (kind=RKIND), dimension(Numcols,Plev), intent(out) :: Uy
+ real (kind=RKIND), dimension(Numcols,Plev), intent(out) :: Omega
+ real (kind=RKIND), dimension(Numcols,Plev,Pcnst), intent(out) :: Tracer
+
+ type (block_type), pointer :: block
+
write(0,*) 'Called MPAS_TO_CAM'
+
+ block => domain % blocklist
+
+ !
+ ! Perform basic sanity check on expected and available field dimensions
+ !
+ if (Numcols /= block % mesh % nCellsSolve) then
+ write(0,*) 'Error: mismatch between Numcols and nCellsSolve: ', Numcols, block % mesh % nCellsSolve
+ end if
+ if (Plev /= block % mesh % nVertLevels) then
+ write(0,*) 'Error: mismatch between Plev and nVertLevels: ', Plev, block % mesh % nVertLevels
+ end if
+ if (Pcnst /= num_scalars) then
+ write(0,*) 'Error: mismatch between Pcnst and num_scalars: ', Pcnst, num_scalars
+ end if
+
+
+ !
+ ! Fill in CAM arrays from block % time_levs(2) % state arrays
+ !
+ Psd = 3.1415926
+ Phis = 3.1415926
+ T = 3.1415926
+ Ux = 3.1415926
+ Uy = 3.1415926
+ Omega = 3.1415926
+ Tracer = 3.1415926
end subroutine mpas_to_cam
- subroutine cam_to_mpas()
+ subroutine cam_to_mpas(Numcols, Plev, Pcnst, T_tend, Ux_tend, Uy_tend, Tracer)
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! SUBROUTINE CAM_TO_MPAS
!
@@ -182,7 +345,42 @@
implicit none
+ integer, intent(in) :: Numcols
+ integer, intent(in) :: Plev
+ integer, intent(in) :: Pcnst
+ real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: T_tend
+ real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: Ux_tend
+ real (kind=RKIND), dimension(Numcols,Plev), intent(in) :: Uy_tend
+ real (kind=RKIND), dimension(Numcols,Plev,Pcnst), intent(in) :: Tracer
+
+ type (block_type), pointer :: block
+
write(0,*) 'Called CAM_TO_MPAS'
+
+ block => domain % blocklist
+
+ !
+ ! Perform basic sanity check on expected and available field dimensions
+ !
+ if (Numcols /= block % mesh % nCellsSolve) then
+ write(0,*) 'Error: mismatch between Numcols and nCellsSolve: ', Numcols, block % mesh % nCellsSolve
+ end if
+ if (Plev /= block % mesh % nVertLevels) then
+ write(0,*) 'Error: mismatch between Plev and nVertLevels: ', Plev, block % mesh % nVertLevels
+ end if
+ if (Pcnst /= num_scalars) then
+ write(0,*) 'Error: mismatch between Pcnst and num_scalars: ', Pcnst, num_scalars
+ end if
+
+
+ !
+ ! Fill in MPAS tendency arrays from arguments
+ !
+
+
+ !
+ ! It might also be necessary to do a ghost-cell update on tendency arrays
+ !
end subroutine cam_to_mpas
Modified: branches/mpas_cam_coupling/src/driver_cam_interface/test_driver.F
===================================================================
--- branches/mpas_cam_coupling/src/driver_cam_interface/test_driver.F        2010-03-17 20:33:00 UTC (rev 144)
+++ branches/mpas_cam_coupling/src/driver_cam_interface/test_driver.F        2010-03-18 17:47:52 UTC (rev 145)
@@ -4,25 +4,57 @@
implicit none
+ integer :: nCellsSolve, maxEdges, nVertLevels, nScalars
+ real (kind=RKIND) :: Pref, Ptop
+ real (kind=RKIND), dimension(:), pointer :: Psd, Phis, Ae, Be, Ac, Bc
+ real (kind=RKIND), dimension(:,:), pointer :: T, Ux, Uy, Omega
+ real (kind=RKIND), dimension(:,:), pointer :: T_tend, Ux_tend, Uy_tend
+ real (kind=RKIND), dimension(:,:,:), pointer :: U, Tracer
integer :: mpi_ierr
call MPI_Init(mpi_ierr)
+ Pref = 100000.
+ Ptop = 21900.
+ nCellsSolve = 100
+ maxEdges = 10
+ nVertLevels = 26
+ nScalars = 3
+
+ allocate(Psd(nCellsSolve))
+ allocate(Phis(nCellsSolve))
+ allocate(Ae(nVertLevels+1))
+ allocate(Be(nVertLevels+1))
+ allocate(Ac(nVertLevels))
+ allocate(Bc(nVertLevels))
+ allocate(T(nCellsSolve,nVertLevels))
+ allocate(U(nCellsSolve,maxEdges,nVertLevels))
+ allocate(Ux(nCellsSolve,nVertLevels))
+ allocate(Uy(nCellsSolve,nVertLevels))
+ allocate(Omega(nCellsSolve,nVertLevels))
+ allocate(Tracer(nCellsSolve,nVertLevels,nScalars))
+ allocate(T_tend(nCellsSolve,nVertLevels))
+ allocate(Ux_tend(nCellsSolve,nVertLevels))
+ allocate(Uy_tend(nCellsSolve,nVertLevels))
+
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! Begin calling CAM-MAS interface routines
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
call mpas_init1(MPI_COMM_WORLD)
- call cam_inidat_to_mpas()
+ call cam_inidat_to_mpas(nCellsSolve, maxEdges, nVertLevels, nScalars, Pref, Ptop, &
+ Ae, Be, Ac, Bc, Psd, Phis, T, U, Tracer)
call mpas_init2()
call mpas_dyn_run()
- call mpas_to_cam()
+ call mpas_to_cam(nCellsSolve, nVertLevels, nScalars, &
+ Psd, Phis, T, Ux, Uy, Omega, Tracer)
- call cam_to_mpas()
+ call cam_to_mpas(nCellsSolve, nVertLevels, nScalars, &
+ T_tend, Ux_tend, Uy_tend, Tracer)
call mpas_final()
@@ -30,6 +62,22 @@
! Finish calling CAM-MAS interface routines
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
+ deallocate(Psd)
+ deallocate(Phis)
+ deallocate(Ae)
+ deallocate(Be)
+ deallocate(Ac)
+ deallocate(Bc)
+ deallocate(T)
+ deallocate(U)
+ deallocate(Ux)
+ deallocate(Uy)
+ deallocate(Omega)
+ deallocate(Tracer)
+ deallocate(T_tend)
+ deallocate(Ux_tend)
+ deallocate(Uy_tend)
+
call MPI_Finalize(mpi_ierr)
stop
</font>
</pre>