Line data Source code
1 : module phys_grid
2 :
3 : !------------------------------------------------------------------------------
4 : !
5 : ! The phys_grid module represents the CAM physics decomposition.
6 : !
7 : ! phys_grid_init receives the physics column info (area, weight, centers)
8 : ! from the dycore.
9 : ! The routine then creates the physics decomposition which
10 : ! is the arrangement of columns across the atmosphere model's
11 : ! MPI tasks as well as the arrangement into groups to
12 : ! facilitate efficient threading.
13 : ! The routine then creates a grid object to allow for data
14 : ! to be read into and written from this decomposition.
15 : ! The phys_grid module also provides interfaces for retrieving information
16 : ! about the decomposition
17 : !
18 : ! Note: This current implementation does not perform load balancing,
19 : ! physics columns ae always on the same task as the corresponding
20 : ! column received from the dycore.
21 : !
22 : !------------------------------------------------------------------------------
23 : use shr_kind_mod, only: r8 => shr_kind_r8
24 : use ppgrid, only: begchunk, endchunk
25 : use physics_column_type, only: physics_column_t
26 : use perf_mod, only: t_adj_detailf, t_startf, t_stopf
27 :
28 : implicit none
29 : private
30 : save
31 :
32 : !!XXgoldyXX: v This needs to be removed to complete the weak scaling transition.
33 : public :: SCATTER_FIELD_TO_CHUNK
34 : !!XXgoldyXX: ^ This needs to be removed to complete the weak scaling transition.
35 :
36 : ! Physics grid management
37 : public :: phys_grid_init ! initialize the physics grid
38 : public :: phys_grid_readnl ! Read the phys_grid_nl namelist
39 : public :: phys_grid_initialized
40 : ! Local task interfaces
41 : public :: get_nlcols_p ! Number of local columns
42 : public :: get_area_p ! area of a physics column in radians squared
43 : public :: get_wght_p ! weight of a physics column in radians squared
44 : public :: get_rlat_p ! latitude of a physics column in radians
45 : public :: get_rlon_p ! longitude of a physics column in radians
46 : public :: get_rlat_all_p ! latitudes of physics cols in chunk (radians)
47 : public :: get_rlon_all_p ! longitudes of physics cols in chunk (radians)
48 : public :: get_lat_p ! latitude of a physics column in degrees
49 : public :: get_lon_p ! longitude of a physics column in degrees
50 : public :: get_lat_all_p ! latitudes of physics cols in chunk (degrees)
51 : public :: get_lon_all_p ! longitudes of physics cols in chunk (degrees)
52 : public :: get_area_all_p ! areas of physics cols in chunk
53 : public :: get_wght_all_p ! weights of physics cols in chunk
54 : public :: get_ncols_p ! number of columns in a chunk
55 : public :: get_gcol_p ! global column index of a physics column
56 : public :: get_gcol_all_p ! global col index of all phys cols in a chunk
57 : public :: get_dyn_col_p ! dynamics local blk number and blk offset(s)
58 : public :: get_chunk_info_p ! chunk index and col # of a physics column
59 : public :: get_grid_dims ! return grid dimensions
60 : ! Physics-dynamics coupling
61 : public :: phys_decomp_to_dyn ! Transfer physics data to dynamics decomp
62 : public :: dyn_decomp_to_phys ! Transfer dynamics data to physics decomp
63 :
64 : ! The identifier for the physics grid
65 : integer, parameter, public :: phys_decomp = 100
66 :
67 : !! PUBLIC TYPES
68 :
69 : ! Physics chunking (thread blocking) data
70 : ! Note that chunks cover local data
71 : type, public :: chunk
72 : integer, private :: ncols = 1 ! # of grid columns in this chunk
73 : integer, private :: chunk_index = -1 ! Local index of this chunk
74 : integer, private, allocatable :: phys_cols(:) ! phys column indices
75 : end type chunk
76 :
77 : !! PRIVATE DATA
78 :
79 : ! dynamics field grid information
80 : ! hdim1_d and hdim2_d are dimensions of rectangular horizontal grid
81 : ! data structure, If 1D data structure, then hdim2_d == 1.
82 : integer :: hdim1_d, hdim2_d
83 :
84 : ! Physics decomposition information
85 : type(physics_column_t), allocatable :: phys_columns(:)
86 :
87 : type(chunk), private, pointer :: chunks(:) => NULL() ! (begchunk:endchunk)
88 :
89 : logical :: phys_grid_set = .false.
90 :
91 : logical :: calc_memory_increase = .false.
92 :
93 : interface get_dyn_col_p
94 : module procedure :: get_dyn_col_p_chunk
95 : module procedure :: get_dyn_col_p_index
96 : end interface get_dyn_col_p
97 :
98 : ! Private interfaces
99 : private :: chunk_info_to_index_p
100 :
101 : !!XXgoldyXX: v temporary interface to allow old code to compile
102 : interface get_lat_all_p
103 : module procedure :: get_lat_all_p_r8 ! The new version
104 : module procedure :: get_lat_all_p_int ! calls endun
105 : end interface get_lat_all_p
106 :
107 : interface get_lon_all_p
108 : module procedure :: get_lon_all_p_r8 ! The new version
109 : module procedure :: get_lon_all_p_int ! calls endun
110 : end interface get_lon_all_p
111 : !!XXgoldyXX: ^ temporary interface to allow old code to compile
112 :
113 :
114 : integer, protected, public :: pver = 0
115 : integer, protected, public :: pverp = 0
116 : integer, protected, public :: num_global_phys_cols = 0
117 : integer, protected, public :: columns_on_task = 0
118 : integer, protected, public :: index_top_layer = 0
119 : integer, protected, public :: index_bottom_layer = 0
120 : integer, protected, public :: index_top_interface = 1
121 : integer, protected, public :: index_bottom_interface = 0
122 :
123 : !==============================================================================
124 : CONTAINS
125 : !==============================================================================
126 :
127 1536 : subroutine phys_grid_readnl(nlfile)
128 : use cam_abortutils, only: endrun
129 : use namelist_utils, only: find_group_name
130 : use cam_logfile, only: iulog
131 : use spmd_utils, only: mpicom, mstrid=>masterprocid, masterproc
132 : use spmd_utils, only: mpi_integer
133 : use ppgrid, only: pcols
134 :
135 : character(len=*), intent(in) :: nlfile
136 :
137 : ! Local variables
138 : integer :: unitn, ierr
139 : character(len=*), parameter :: sub = 'phys_grid_readnl'
140 :
141 : integer :: phys_alltoall = -HUGE(1)
142 : integer :: phys_loadbalance = -HUGE(1)
143 : integer :: phys_twin_algorithm = -HUGE(1)
144 : integer :: phys_chnk_per_thd = -HUGE(1)
145 :
146 : namelist /phys_grid_nl/ phys_alltoall, phys_loadbalance, &
147 : phys_twin_algorithm, phys_chnk_per_thd
148 : !------------------------------------------------------------------------
149 :
150 : ! Read namelist
151 1536 : if (masterproc) then
152 2 : open(newunit=unitn, file=trim(nlfile), status='old')
153 2 : call find_group_name(unitn, 'phys_grid_nl', status=ierr)
154 2 : if (ierr == 0) then
155 0 : read(unitn, phys_grid_nl, iostat=ierr)
156 0 : if (ierr /= 0) then
157 0 : call endrun(sub//': FATAL: reading namelist')
158 : end if
159 : end if
160 2 : close(unitn)
161 : end if
162 :
163 1536 : call mpi_bcast(phys_alltoall, 1, mpi_integer, mstrid, mpicom, ierr)
164 1536 : call mpi_bcast(phys_loadbalance, 1, mpi_integer, mstrid, mpicom, ierr)
165 1536 : call mpi_bcast(phys_twin_algorithm, 1, mpi_integer, mstrid, mpicom, ierr)
166 1536 : call mpi_bcast(phys_chnk_per_thd, 1, mpi_integer, mstrid, mpicom, ierr)
167 :
168 1536 : if (masterproc) then
169 2 : write(iulog,*) 'PHYS_GRID options:'
170 2 : write(iulog,*) ' Using PCOLS =', pcols
171 2 : write(iulog,*) ' phys_loadbalance = (not used)'
172 2 : write(iulog,*) ' phys_twin_algorithm = (not used)'
173 2 : write(iulog,*) ' phys_alltoall = (not used)'
174 2 : write(iulog,*) ' chunks_per_thread = (not used)'
175 : end if
176 :
177 1536 : end subroutine phys_grid_readnl
178 :
179 : !========================================================================
180 :
181 1536 : subroutine phys_grid_init()
182 : use mpi, only: MPI_INTEGER, MPI_REAL8, MPI_MIN, MPI_MAX
183 : use shr_mem_mod, only: shr_mem_getusage
184 : use cam_abortutils, only: endrun
185 : use cam_logfile, only: iulog
186 : use spmd_utils, only: npes, mpicom, masterprocid, masterproc, iam
187 : use ppgrid, only: pcols
188 : use dyn_grid, only: get_dyn_grid_info, physgrid_copy_attributes_d
189 : use cam_grid_support, only: cam_grid_register, cam_grid_attribute_register
190 : use cam_grid_support, only: iMap, hclen => max_hcoordname_len
191 : use cam_grid_support, only: horiz_coord_t, horiz_coord_create
192 : use cam_grid_support, only: cam_grid_attribute_copy, cam_grid_attr_exists
193 : use shr_const_mod, only: PI => SHR_CONST_PI
194 :
195 : ! Local variables
196 : integer :: index
197 : integer :: col_index, phys_col
198 : integer :: ichnk, icol, ncol, gcol
199 : integer :: num_chunks
200 1536 : type(physics_column_t), allocatable :: dyn_columns(:) ! Dyn decomp
201 : ! Maps and values for physics grid
202 1536 : real(r8), pointer :: lonvals(:)
203 1536 : real(r8), pointer :: latvals(:)
204 : real(r8) :: lonmin, latmin
205 1536 : integer(iMap), pointer :: grid_map(:,:)
206 1536 : integer(iMap), allocatable :: coord_map(:)
207 : type(horiz_coord_t), pointer :: lat_coord
208 : type(horiz_coord_t), pointer :: lon_coord
209 1536 : real(r8), pointer :: area_d(:)
210 1536 : real(r8), pointer :: areawt_d(:)
211 : real(r8) :: mem_hw_beg, mem_hw_end
212 : real(r8) :: mem_beg, mem_end
213 : logical :: unstructured
214 : real(r8) :: temp ! For MPI
215 : integer :: ierr ! For MPI
216 1536 : character(len=hclen), pointer :: copy_attributes(:)
217 : character(len=hclen) :: copy_gridname
218 : character(len=*), parameter :: subname = 'phys_grid_init: '
219 : real(r8), parameter :: rarea_sphere = 1.0_r8 / (4.0_r8*PI)
220 :
221 1536 : nullify(lonvals)
222 1536 : nullify(latvals)
223 1536 : nullify(grid_map)
224 1536 : nullify(lat_coord)
225 1536 : nullify(lon_coord)
226 1536 : nullify(area_d)
227 1536 : nullify(areawt_d)
228 1536 : nullify(copy_attributes)
229 :
230 1536 : if (calc_memory_increase) then
231 0 : call shr_mem_getusage(mem_hw_beg, mem_beg)
232 : end if
233 :
234 1536 : call t_adj_detailf(-2)
235 1536 : call t_startf("phys_grid_init")
236 :
237 : ! Gather info from the dycore
238 : call get_dyn_grid_info(hdim1_d, hdim2_d, pver, index_top_layer, &
239 1536 : index_bottom_layer, unstructured, dyn_columns)
240 : ! hdim1_d * hdim2_d is the total number of columns
241 1536 : num_global_phys_cols = hdim1_d * hdim2_d
242 1536 : pverp = pver + 1
243 : !!XXgoldyXX: Can we enforce interface numbering separate from dycore?
244 : !!XXgoldyXX: This will work for both CAM and WRF/MPAS physics
245 : !!XXgoldyXX: This only has a 50% chance of working on a single level model
246 1536 : if (index_top_layer < index_bottom_layer) then
247 1536 : index_top_interface = index_top_layer
248 1536 : index_bottom_interface = index_bottom_layer + 1
249 : else
250 0 : index_bottom_interface = index_bottom_layer
251 0 : index_top_interface = index_top_layer + 1
252 : end if
253 :
254 : ! Set up the physics decomposition
255 1536 : columns_on_task = size(dyn_columns)
256 1536 : if (allocated(phys_columns)) then
257 0 : deallocate(phys_columns)
258 : end if
259 104880 : allocate(phys_columns(columns_on_task))
260 1536 : if (columns_on_task > 0) then
261 1536 : col_index = columns_on_task
262 1536 : num_chunks = col_index / pcols
263 1536 : if ((num_chunks * pcols) < col_index) then
264 1536 : num_chunks = num_chunks + 1
265 : end if
266 1536 : begchunk = 1
267 1536 : endchunk = begchunk + num_chunks - 1
268 : else
269 : ! We do not support tasks with no physics columns
270 0 : call endrun(subname//'No columns on task, use fewer tasks')
271 : end if
272 10800 : allocate(chunks(begchunk:endchunk))
273 1536 : col_index = 0
274 : ! Simple chunk assignment
275 7728 : do index = begchunk, endchunk
276 6192 : chunks(index)%ncols = MIN(pcols, (columns_on_task - col_index))
277 6192 : chunks(index)%chunk_index = index
278 18576 : allocate(chunks(index)%phys_cols(chunks(index)%ncols))
279 104928 : do phys_col = 1, chunks(index)%ncols
280 97200 : col_index = col_index + 1
281 : ! Copy information supplied by the dycore
282 97200 : phys_columns(col_index) = dyn_columns(col_index)
283 : ! Fill in physics decomp info
284 97200 : phys_columns(col_index)%phys_task = iam
285 97200 : phys_columns(col_index)%local_phys_chunk = index
286 97200 : phys_columns(col_index)%phys_chunk_index = phys_col
287 103392 : chunks(index)%phys_cols(phys_col) = col_index
288 : end do
289 : end do
290 :
291 98736 : deallocate(dyn_columns)
292 :
293 : ! Add physics-package grid to set of CAM grids
294 : ! physgrid always uses 'lat' and 'lon' as coordinate names; If dynamics
295 : ! grid is different, it will use different coordinate names
296 :
297 : ! First, create a map for the physics grid
298 : ! It's structure will depend on whether or not the physics grid is
299 : ! unstructured
300 1536 : if (unstructured) then
301 4608 : allocate(grid_map(3, pcols * (endchunk - begchunk + 1)))
302 : else
303 0 : allocate(grid_map(4, pcols * (endchunk - begchunk + 1)))
304 : end if
305 397824 : grid_map = 0_iMap
306 4608 : allocate(latvals(size(grid_map, 2)))
307 3072 : allocate(lonvals(size(grid_map, 2)))
308 :
309 1536 : lonmin = 1000.0_r8 ! Out of longitude range
310 1536 : latmin = 1000.0_r8 ! Out of latitude range
311 1536 : index = 0
312 7728 : do ichnk = begchunk, endchunk
313 6192 : ncol = chunks(ichnk)%ncols ! Too soon to call get_ncols_p
314 106800 : do icol = 1, pcols
315 99072 : index = index + 1
316 99072 : if (icol <= ncol) then
317 97200 : col_index = chunks(ichnk)%phys_cols(icol)
318 97200 : latvals(index) = phys_columns(col_index)%lat_deg
319 97200 : if (latvals(index) < latmin) then
320 9410 : latmin = latvals(index)
321 : end if
322 97200 : lonvals(index) = phys_columns(col_index)%lon_deg
323 97200 : if (lonvals(index) < lonmin) then
324 6592 : lonmin = lonvals(index)
325 : end if
326 : else
327 1872 : col_index = -1
328 1872 : latvals(index) = 1000.0_r8
329 1872 : lonvals(index) = 1000.0_r8
330 : end if
331 99072 : grid_map(1, index) = int(icol, iMap)
332 99072 : grid_map(2, index) = int(ichnk, iMap)
333 105264 : if (icol <= ncol) then
334 97200 : if (unstructured) then
335 97200 : gcol = phys_columns(col_index)%global_col_num
336 97200 : if (gcol > 0) then
337 97200 : grid_map(3, index) = int(gcol, iMap)
338 : end if ! else entry remains 0
339 : else
340 : ! lon
341 0 : gcol = phys_columns(col_index)%coord_indices(1)
342 0 : if (gcol > 0) then
343 0 : grid_map(3, index) = int(gcol, iMap)
344 : end if ! else entry remains 0
345 : ! lat
346 0 : gcol = phys_columns(col_index)%coord_indices(2)
347 0 : if (gcol > 0) then
348 0 : grid_map(4, index) = gcol
349 : end if ! else entry remains 0
350 : end if
351 : end if ! Else entry remains 0
352 : end do
353 : end do
354 :
355 : ! Note that if the dycore is using the same points as the physics grid,
356 : ! it will have already set up 'lat' and 'lon' axes for
357 : ! the physics grid
358 : ! However, these will be in the dynamics decomposition
359 :
360 1536 : if (unstructured) then
361 : lon_coord => horiz_coord_create('lon', 'ncol', num_global_phys_cols, &
362 : 'longitude', 'degrees_east', 1, size(lonvals), lonvals, &
363 1536 : map=grid_map(3,:))
364 : lat_coord => horiz_coord_create('lat', 'ncol', num_global_phys_cols, &
365 : 'latitude', 'degrees_north', 1, size(latvals), latvals, &
366 1536 : map=grid_map(3,:))
367 : else
368 0 : allocate(coord_map(size(grid_map, 2)))
369 : ! We need a global minimum longitude and latitude
370 0 : if (npes > 1) then
371 0 : temp = lonmin
372 : call MPI_allreduce(temp, lonmin, 1, MPI_INTEGER, MPI_MIN, &
373 0 : mpicom, ierr)
374 0 : temp = latmin
375 : call MPI_allreduce(temp, latmin, 1, MPI_INTEGER, MPI_MIN, &
376 0 : mpicom, ierr)
377 : ! Create lon coord map which only writes from one of each unique lon
378 0 : where(latvals == latmin)
379 0 : coord_map(:) = grid_map(3, :)
380 : elsewhere
381 : coord_map(:) = 0_iMap
382 : end where
383 : lon_coord => horiz_coord_create('lon', 'lon', hdim1_d, &
384 : 'longitude', 'degrees_east', 1, size(lonvals), lonvals, &
385 0 : map=coord_map)
386 :
387 : ! Create lat coord map which only writes from one of each unique lat
388 0 : where(lonvals == lonmin)
389 0 : coord_map(:) = grid_map(4, :)
390 : elsewhere
391 : coord_map(:) = 0_iMap
392 : end where
393 : lat_coord => horiz_coord_create('lat', 'lat', hdim2_d, &
394 : 'latitude', 'degrees_north', 1, size(latvals), latvals, &
395 0 : map=coord_map)
396 0 : deallocate(coord_map)
397 : end if
398 : end if
399 : call cam_grid_register('physgrid', phys_decomp, lat_coord, lon_coord, &
400 1536 : grid_map, unstruct=unstructured, block_indexed=.true.)
401 : ! Copy required attributes from the dynamics array
402 1536 : nullify(copy_attributes)
403 1536 : call physgrid_copy_attributes_d(copy_gridname, copy_attributes)
404 4608 : do index = 1, size(copy_attributes)
405 : call cam_grid_attribute_copy(copy_gridname, 'physgrid', &
406 4608 : copy_attributes(index))
407 : end do
408 :
409 1536 : if (.not. cam_grid_attr_exists('physgrid', 'area')) then
410 : ! Physgrid always needs an area attribute.
411 1536 : if (unstructured) then
412 : ! If we did not inherit one from the dycore (i.e., physics and
413 : ! dynamics are on different grids), create that attribute here
414 : ! (Note, a separate physics grid is only supported for
415 : ! unstructured grids).
416 4608 : allocate(area_d(size(grid_map, 2)))
417 98736 : do col_index = 1, columns_on_task
418 98736 : area_d(col_index) = phys_columns(col_index)%area
419 : end do
420 : call cam_grid_attribute_register('physgrid', 'area', &
421 1536 : 'physics column areas', 'ncol', area_d, map=grid_map(3,:))
422 1536 : nullify(area_d) ! Belongs to attribute now
423 :
424 4608 : allocate(areawt_d(size(grid_map, 2)))
425 98736 : do col_index = 1, columns_on_task
426 98736 : areawt_d(col_index) = phys_columns(col_index)%weight*rarea_sphere
427 : end do
428 : call cam_grid_attribute_register('physgrid', 'areawt', &
429 1536 : 'physics column area weight', 'ncol', areawt_d, map=grid_map(3,:))
430 1536 : nullify(areawt_d) ! Belongs to attribute now
431 : else
432 0 : call endrun(subname//"No 'area' attribute from dycore")
433 : end if
434 : end if
435 : ! Cleanup pointers (they belong to the grid now)
436 1536 : nullify(grid_map)
437 1536 : deallocate(latvals)
438 : nullify(latvals)
439 1536 : deallocate(lonvals)
440 : nullify(lonvals)
441 : ! Cleanup, we are responsible for copy attributes
442 1536 : if (associated(copy_attributes)) then
443 1536 : deallocate(copy_attributes)
444 : nullify(copy_attributes)
445 : end if
446 :
447 : ! Set flag indicating physics grid is now set
448 1536 : phys_grid_set = .true.
449 :
450 1536 : call t_stopf("phys_grid_init")
451 1536 : call t_adj_detailf(+2)
452 :
453 1536 : if (calc_memory_increase) then
454 0 : call shr_mem_getusage(mem_hw_end, mem_end)
455 0 : temp = mem_end - mem_beg
456 : call MPI_reduce(temp, mem_end, 1, MPI_REAL8, MPI_MAX, masterprocid, &
457 0 : mpicom, ierr)
458 0 : if (masterproc) then
459 0 : write(iulog, *) 'phys_grid_init: Increase in memory usage = ', &
460 0 : mem_end, ' (MB)'
461 : end if
462 0 : temp = mem_hw_end - mem_hw_beg
463 : call MPI_reduce(temp, mem_hw_end, 1, MPI_REAL8, MPI_MAX, &
464 0 : masterprocid, mpicom, ierr)
465 0 : if (masterproc) then
466 0 : write(iulog, *) subname, 'Increase in memory highwater = ', &
467 0 : mem_end, ' (MB)'
468 : end if
469 : end if
470 :
471 3072 : end subroutine phys_grid_init
472 :
473 : !========================================================================
474 :
475 2964600 : integer function chunk_info_to_index_p(lcid, col, subname_in)
476 1536 : use cam_logfile, only: iulog
477 : use cam_abortutils, only: endrun
478 : ! Return the physics column index indicated by
479 : ! <lcid> (chunk) and <col> (column).
480 :
481 : ! Dummy arguments
482 : integer, intent(in) :: lcid ! local chunk id
483 : integer, intent(in) :: col ! Column index
484 : character(len=*), optional, intent(in) :: subname_in
485 : ! Local variables
486 : character(len=128) :: errmsg
487 : character(len=*), parameter :: subname = 'chunk_info_to_index_p: '
488 :
489 2964600 : if (.not. phys_grid_initialized()) then
490 0 : if (present(subname_in)) then
491 0 : call endrun(trim(subname_in)//'physics grid not initialized')
492 : else
493 0 : call endrun(subname//'physics grid not initialized')
494 : end if
495 2964600 : else if ((lcid < begchunk) .or. (lcid > endchunk)) then
496 0 : if (present(subname_in)) then
497 0 : write(errmsg, '(a,3(a,i0))') trim(subname_in), 'lcid (', lcid, &
498 0 : ') out of range (', begchunk, ' to ', endchunk
499 : else
500 0 : write(errmsg, '(a,3(a,i0))') subname, 'lcid (', lcid, &
501 0 : ') out of range (', begchunk, ' to ', endchunk
502 : end if
503 0 : write(iulog, *) trim(errmsg)
504 0 : call endrun(trim(errmsg))
505 2964600 : else if ((col < 1) .or. (col > get_ncols_p(lcid))) then
506 0 : if (present(subname_in)) then
507 0 : write(errmsg, '(a,2(a,i0))') trim(subname_in), 'col (', col, &
508 0 : ') out of range (1 to ', get_ncols_p(lcid)
509 : else
510 0 : write(errmsg, '(a,2(a,i0))') subname, 'col (', col, &
511 0 : ') out of range (1 to ', get_ncols_p(lcid)
512 : end if
513 0 : write(iulog, *) trim(errmsg)
514 0 : call endrun(trim(errmsg))
515 : end if
516 2964600 : chunk_info_to_index_p = chunks(lcid)%phys_cols(col)
517 2964600 : end function chunk_info_to_index_p
518 :
519 : !========================================================================
520 :
521 16078824 : logical function phys_grid_initialized()
522 : ! Return .true. if the physics grid is initialized, otherwise .false.
523 16078824 : phys_grid_initialized = phys_grid_set
524 16078824 : end function phys_grid_initialized
525 :
526 : !========================================================================
527 :
528 104448 : integer function get_nlcols_p()
529 104448 : get_nlcols_p = columns_on_task
530 104448 : end function get_nlcols_p
531 :
532 : !========================================================================
533 :
534 0 : real(r8) function get_rlat_p(lcid, col)
535 : !-----------------------------------------------------------------------
536 : !
537 : ! get_rlat_p: latitude of a physics column in radians
538 : !
539 : !-----------------------------------------------------------------------
540 :
541 : ! Dummy argument
542 : integer, intent(in) :: lcid
543 : integer, intent(in) :: col
544 : ! Local variables
545 : integer :: index
546 : character(len=*), parameter :: subname = 'get_rlat_p'
547 :
548 0 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
549 0 : get_rlat_p = phys_columns(index)%lat_rad
550 :
551 0 : end function get_rlat_p
552 :
553 : !========================================================================
554 :
555 0 : real(r8) function get_rlon_p(lcid, col)
556 : !-----------------------------------------------------------------------
557 : !
558 : ! get_rlon_p: longitude of a physics column in radians
559 : !
560 : !-----------------------------------------------------------------------
561 :
562 : ! Dummy argument
563 : integer, intent(in) :: lcid
564 : integer, intent(in) :: col
565 : ! Local variables
566 : integer :: index
567 : character(len=*), parameter :: subname = 'get_rlon_p'
568 :
569 0 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
570 0 : get_rlon_p = phys_columns(index)%lon_rad
571 :
572 0 : end function get_rlon_p
573 :
574 : !========================================================================
575 :
576 774000 : subroutine get_rlat_all_p(lcid, rlatdim, rlats)
577 : use cam_abortutils, only: endrun
578 : !-----------------------------------------------------------------------
579 : !
580 : ! get_rlat_all_p: Return all latitudes (in radians) for chunk, <lcid>
581 : !
582 : !-----------------------------------------------------------------------
583 : ! Dummy Arguments
584 : integer, intent(in) :: lcid ! local chunk id
585 : integer, intent(in) :: rlatdim ! declared size of output array
586 : real(r8), intent(out) :: rlats(rlatdim) ! array of latitudes
587 :
588 : ! Local variables
589 : integer :: index ! loop index
590 : integer :: phys_ind
591 : character(len=*), parameter :: subname = 'get_rlat_all_p: '
592 :
593 : !-----------------------------------------------------------------------
594 774000 : if ((lcid < begchunk) .or. (lcid > endchunk)) then
595 0 : call endrun(subname//'chunk index out of range')
596 : end if
597 12924000 : do index = 1, MIN(get_ncols_p(lcid, subname_in=subname), rlatdim)
598 12150000 : phys_ind = chunks(lcid)%phys_cols(index)
599 12924000 : rlats(index) = phys_columns(phys_ind)%lat_rad
600 : end do
601 :
602 774000 : end subroutine get_rlat_all_p
603 :
604 : !========================================================================
605 :
606 585144 : subroutine get_rlon_all_p(lcid, rlondim, rlons)
607 : use cam_abortutils, only: endrun
608 : !-----------------------------------------------------------------------
609 : !
610 : ! get_rlon_all_p:: Return all longitudes (in radians) for chunk, <lcid>
611 : !
612 : !-----------------------------------------------------------------------
613 : ! Dummy Arguments
614 : integer, intent(in) :: lcid ! local chunk id
615 : integer, intent(in) :: rlondim ! declared size of output array
616 : real(r8), intent(out) :: rlons(rlondim) ! array of longitudes
617 :
618 : ! Local variables
619 : integer :: index ! loop index
620 : integer :: phys_ind
621 : character(len=*), parameter :: subname = 'get_rlon_all_p: '
622 :
623 : !-----------------------------------------------------------------------
624 585144 : if ((lcid < begchunk) .or. (lcid > endchunk)) then
625 0 : call endrun(subname//'chunk index out of range')
626 : end if
627 9770544 : do index = 1, MIN(get_ncols_p(lcid, subname_in=subname), rlondim)
628 9185400 : phys_ind = chunks(lcid)%phys_cols(index)
629 9770544 : rlons(index) = phys_columns(phys_ind)%lon_rad
630 : end do
631 :
632 585144 : end subroutine get_rlon_all_p
633 :
634 : !========================================================================
635 :
636 0 : real(r8) function get_lat_p(lcid, col)
637 : !-----------------------------------------------------------------------
638 : !
639 : ! get_lat_p: latitude of a physics column in degrees
640 : !
641 : !-----------------------------------------------------------------------
642 :
643 : ! Dummy argument
644 : integer, intent(in) :: lcid
645 : integer, intent(in) :: col
646 : ! Local variables
647 : integer :: index
648 : character(len=*), parameter :: subname = 'get_lat_p'
649 :
650 0 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
651 0 : get_lat_p = phys_columns(index)%lat_deg
652 :
653 0 : end function get_lat_p
654 :
655 : !========================================================================
656 :
657 0 : real(r8) function get_lon_p(lcid, col)
658 : !-----------------------------------------------------------------------
659 : !
660 : ! get_lon_p: longitude of a physics column in degrees
661 : !
662 : !-----------------------------------------------------------------------
663 :
664 : ! Dummy argument
665 : integer, intent(in) :: lcid
666 : integer, intent(in) :: col
667 : ! Local variables
668 : integer :: index
669 : character(len=*), parameter :: subname = 'get_lon_p'
670 :
671 0 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
672 0 : get_lon_p = phys_columns(index)%lon_deg
673 :
674 0 : end function get_lon_p
675 :
676 : !========================================================================
677 :
678 0 : subroutine get_lat_all_p_r8(lcid, latdim, lats)
679 : use cam_abortutils, only: endrun
680 : !-----------------------------------------------------------------------
681 : !
682 : ! get_lat_all_p: Return all latitudes (in degrees) for chunk, <lcid>
683 : !
684 : !-----------------------------------------------------------------------
685 : ! Dummy Arguments
686 : integer, intent(in) :: lcid ! local chunk id
687 : integer, intent(in) :: latdim ! declared size of output array
688 : real(r8), intent(out) :: lats(latdim) ! array of latitudes
689 :
690 : ! Local variables
691 : integer :: index ! loop index
692 : integer :: phys_ind
693 : character(len=*), parameter :: subname = 'get_lat_all_p: '
694 :
695 : !-----------------------------------------------------------------------
696 0 : if ((lcid < begchunk) .or. (lcid > endchunk)) then
697 0 : call endrun(subname//'chunk index out of range')
698 : end if
699 0 : do index = 1, MIN(get_ncols_p(lcid, subname_in=subname), latdim)
700 0 : phys_ind = chunks(lcid)%phys_cols(index)
701 0 : lats(index) = phys_columns(phys_ind)%lat_deg
702 : end do
703 :
704 0 : end subroutine get_lat_all_p_r8
705 :
706 : !========================================================================
707 :
708 0 : subroutine get_lon_all_p_r8(lcid, londim, lons)
709 : use cam_abortutils, only: endrun
710 : !-----------------------------------------------------------------------
711 : !
712 : ! get_lon_all_p:: Return all longitudes (in degrees) for chunk, <lcid>
713 : !
714 : !-----------------------------------------------------------------------
715 : ! Dummy Arguments
716 : integer, intent(in) :: lcid ! local chunk id
717 : integer, intent(in) :: londim ! declared size of output array
718 : real(r8), intent(out) :: lons(londim) ! array of longitudes
719 :
720 : ! Local variables
721 : integer :: index ! loop index
722 : integer :: phys_ind
723 : character(len=*), parameter :: subname = 'get_lon_all_p: '
724 :
725 : !-----------------------------------------------------------------------
726 0 : if ((lcid < begchunk) .or. (lcid > endchunk)) then
727 0 : call endrun(subname//'chunk index out of range')
728 : end if
729 0 : do index = 1, MIN(get_ncols_p(lcid, subname_in=subname), londim)
730 0 : phys_ind = chunks(lcid)%phys_cols(index)
731 0 : lons(index) = phys_columns(phys_ind)%lon_deg
732 : end do
733 :
734 0 : end subroutine get_lon_all_p_r8
735 :
736 : !========================================================================
737 :
738 117648 : subroutine get_area_all_p(lcid, areadim, areas)
739 : use cam_abortutils, only: endrun
740 : !-----------------------------------------------------------------------
741 : !
742 : ! get_area_all_p: Return all areas for chunk, <lcid>
743 : !
744 : !-----------------------------------------------------------------------
745 : ! Dummy Arguments
746 : integer, intent(in) :: lcid ! local chunk id
747 : integer, intent(in) :: areadim ! declared size of output array
748 : real(r8), intent(out) :: areas(areadim) ! array of areas
749 :
750 : ! Local variables
751 : integer :: index ! loop index
752 : integer :: phys_ind
753 : character(len=*), parameter :: subname = 'get_area_all_p: '
754 :
755 : !-----------------------------------------------------------------------
756 117648 : if ((lcid < begchunk) .or. (lcid > endchunk)) then
757 0 : call endrun(subname//'chunk index out of range')
758 : end if
759 1964448 : do index = 1, MIN(get_ncols_p(lcid, subname_in=subname), areadim)
760 1846800 : phys_ind = chunks(lcid)%phys_cols(index)
761 1964448 : areas(index) = phys_columns(phys_ind)%area
762 : end do
763 :
764 117648 : end subroutine get_area_all_p
765 :
766 : !========================================================================
767 :
768 733752 : subroutine get_wght_all_p(lcid, wghtdim, wghts)
769 : use cam_abortutils, only: endrun
770 : !-----------------------------------------------------------------------
771 : !
772 : ! get_wght_all_p: Return all weights for chunk, <lcid>
773 : !
774 : !-----------------------------------------------------------------------
775 : ! Dummy Arguments
776 : integer, intent(in) :: lcid ! local chunk id
777 : integer, intent(in) :: wghtdim ! declared size of output array
778 : real(r8), intent(out) :: wghts(wghtdim) ! array of weights
779 :
780 : ! Local variables
781 : integer :: index ! loop index
782 : integer :: phys_ind
783 : character(len=*), parameter :: subname = 'get_wght_all_p: '
784 :
785 : !-----------------------------------------------------------------------
786 733752 : if ((lcid < begchunk) .or. (lcid > endchunk)) then
787 0 : call endrun(subname//'chunk index out of range')
788 : end if
789 12251952 : do index = 1, MIN(get_ncols_p(lcid, subname_in=subname), wghtdim)
790 11518200 : phys_ind = chunks(lcid)%phys_cols(index)
791 12251952 : wghts(index) = phys_columns(phys_ind)%weight
792 : end do
793 :
794 733752 : end subroutine get_wght_all_p
795 :
796 : !========================================================================
797 :
798 9215424 : integer function get_ncols_p(lcid, subname_in)
799 : use cam_abortutils, only: endrun
800 : !-----------------------------------------------------------------------
801 : !
802 : ! get_ncols_p: Return number of columns in chunk given the local chunk id.
803 : !
804 : !-----------------------------------------------------------------------
805 : ! Dummy arguments
806 : integer, intent(in) :: lcid ! local chunk id
807 : character(len=*), optional, intent(in) :: subname_in
808 :
809 9215424 : if (.not. phys_grid_initialized()) then
810 0 : if (present(subname_in)) then
811 0 : call endrun(trim(subname_in)//'physics grid not initialized')
812 : else
813 0 : call endrun('get_ncols_p: physics grid not initialized')
814 : end if
815 : else
816 9215424 : get_ncols_p = chunks(lcid)%ncols
817 : end if
818 :
819 9215424 : end function get_ncols_p
820 :
821 : !========================================================================
822 :
823 2770200 : real(r8) function get_area_p(lcid, col)
824 : ! area of a physics column in radians squared
825 :
826 : ! Dummy arguments
827 : integer, intent(in) :: lcid ! Chunk number
828 : integer, intent(in) :: col ! <lcid> column
829 : ! Local variables
830 : integer :: index
831 : character(len=*), parameter :: subname = 'get_area_p'
832 :
833 2770200 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
834 2770200 : get_area_p = phys_columns(index)%area
835 :
836 2770200 : end function get_area_p
837 :
838 : !========================================================================
839 :
840 0 : real(r8) function get_wght_p(lcid, col)
841 : ! weight of a physics column in radians squared
842 :
843 : ! Dummy arguments
844 : integer, intent(in) :: lcid ! Chunk number
845 : integer, intent(in) :: col ! <lcid> column
846 : ! Local variables
847 : integer :: index
848 : character(len=*), parameter :: subname = 'get_wght_p'
849 :
850 0 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
851 0 : get_wght_p = phys_columns(index)%weight
852 :
853 0 : end function get_wght_p
854 :
855 : !========================================================================
856 :
857 194400 : integer function get_gcol_p(lcid, col)
858 : ! global column index of a physics column
859 :
860 : ! Dummy arguments
861 : integer, intent(in) :: lcid ! local chunk id
862 : integer, intent(in) :: col ! column index
863 : ! Local variables
864 : integer :: index
865 : character(len=*), parameter :: subname = 'get_gcol_p: '
866 :
867 194400 : index = chunk_info_to_index_p(lcid, col, subname_in=subname)
868 194400 : get_gcol_p = phys_columns(index)%global_col_num
869 :
870 194400 : end function get_gcol_p
871 :
872 : !========================================================================
873 :
874 0 : subroutine get_dyn_col_p_chunk(lcid, col, blk_num, blk_ind, caller)
875 : use cam_abortutils, only: endrun
876 : ! Return the dynamics local block number and block offset(s) for
877 : ! the physics column indicated by <lcid> (chunk) and <col> (column).
878 :
879 : ! Dummy arguments
880 : integer, intent(in) :: lcid ! local chunk id
881 : integer, intent(in) :: col ! Column index
882 : integer, intent(out) :: blk_num ! Local dynamics block index
883 : integer, intent(out) :: blk_ind(:) ! Local dynamics block offset(s)
884 : character(len=*), optional, intent(in) :: caller ! Calling routine
885 : ! Local variables
886 : integer :: index
887 : integer :: off_size
888 : character(len=*), parameter :: subname = 'get_dyn_col_p_chunk: '
889 :
890 0 : index = chunk_info_to_index_p(lcid, col)
891 0 : off_size = SIZE(phys_columns(index)%dyn_block_index, 1)
892 0 : if (SIZE(blk_ind, 1) < off_size) then
893 0 : if (present(caller)) then
894 0 : call endrun(trim(caller)//': blk_ind too small')
895 : else
896 0 : call endrun(subname//'blk_ind too small')
897 : end if
898 : end if
899 0 : blk_num = phys_columns(index)%local_dyn_block
900 0 : blk_ind(1:off_size) = phys_columns(index)%dyn_block_index(1:off_size)
901 0 : if (SIZE(blk_ind, 1) > off_size) then
902 0 : blk_ind(off_size+1:) = -1
903 : end if
904 :
905 0 : end subroutine get_dyn_col_p_chunk
906 :
907 : !========================================================================
908 :
909 1944000 : subroutine get_dyn_col_p_index(index, blk_num, blk_ind)
910 : use cam_logfile, only: iulog
911 : use cam_abortutils, only: endrun
912 : ! Return the dynamics local block number and block offset(s) for
913 : ! the physics column indicated by <index>.
914 :
915 : ! Dummy arguments
916 : integer, intent(in) :: index ! index of local physics column
917 : integer, intent(out) :: blk_num ! Local dynamics block index
918 : integer, intent(out) :: blk_ind(:) ! Local dynamics block offset(s)
919 : ! Local variables
920 : integer :: off_size
921 : character(len=128) :: errmsg
922 : character(len=*), parameter :: subname = 'get_dyn_col_p_index: '
923 :
924 1944000 : if (.not. phys_grid_initialized()) then
925 0 : call endrun(subname//'physics grid not initialized')
926 1944000 : else if ((index < 1) .or. (index > columns_on_task)) then
927 0 : write(errmsg, '(a,2(a,i0))') subname, 'index (', index, &
928 0 : ') out of range (1 to ', columns_on_task
929 0 : write(iulog, *) trim(errmsg)
930 0 : call endrun(trim(errmsg))
931 : else
932 1944000 : off_size = SIZE(phys_columns(index)%dyn_block_index, 1)
933 1944000 : if (SIZE(blk_ind, 1) < off_size) then
934 0 : call endrun(subname//'blk_ind too small')
935 : end if
936 1944000 : blk_num = phys_columns(index)%local_dyn_block
937 3888000 : blk_ind(1:off_size) = phys_columns(index)%dyn_block_index(1:off_size)
938 1944000 : if (SIZE(blk_ind, 1) > off_size) then
939 0 : blk_ind(off_size+1:) = -1
940 : end if
941 : end if
942 :
943 1944000 : end subroutine get_dyn_col_p_index
944 :
945 : !========================================================================
946 :
947 6192 : subroutine get_gcol_all_p(lcid, gdim, gcols)
948 : use cam_logfile, only: iulog
949 : use cam_abortutils, only: endrun
950 : use spmd_utils, only: masterproc
951 : ! collect global column indices of all physics columns in a chunk
952 :
953 : ! Dummy arguments
954 : integer, intent(in) :: lcid ! local chunk id
955 : integer, intent(in) :: gdim ! gcols dimension
956 : integer, intent(out) :: gcols(:) ! global column indices
957 : ! Local variables
958 : integer :: ncol, col_ind
959 : character(len=128) :: errmsg
960 : character(len=*), parameter :: subname = 'get_gcol_all_p: '
961 :
962 6192 : if (.not. phys_grid_initialized()) then
963 0 : call endrun(subname//'physics grid not initialized')
964 6192 : else if ((lcid < begchunk) .or. (lcid > endchunk)) then
965 0 : write(errmsg, '(a,3(a,i0))') subname, 'lcid (', lcid, &
966 0 : ') out of range (', begchunk, ' to ', endchunk
967 0 : write(iulog, *) trim(errmsg)
968 0 : call endrun(trim(errmsg))
969 : else
970 6192 : ncol = chunks(lcid)%ncols
971 6192 : if (gdim < ncol) then
972 0 : if (masterproc) then
973 0 : write(iulog, '(2a,2(i0,a))') subname, 'WARNING: gdim (', gdim, &
974 0 : ') < ncol (', ncol,'), not all indices will be filled.'
975 : end if
976 0 : gcols(gdim+1:ncol) = -1
977 : end if
978 103392 : do col_ind = 1, MIN(ncol, gdim)
979 103392 : gcols(col_ind) = get_gcol_p(lcid, col_ind)
980 : end do
981 : end if
982 :
983 6192 : end subroutine get_gcol_all_p
984 :
985 : !========================================================================
986 :
987 1944000 : subroutine get_chunk_info_p(index, lchnk, icol)
988 : use cam_logfile, only: iulog
989 : use cam_abortutils, only: endrun
990 : ! local chunk index and column number of a physics column
991 :
992 : ! Dummy arguments
993 : integer, intent(in) :: index
994 : integer, intent(out) :: lchnk
995 : integer, intent(out) :: icol
996 : ! Local variables
997 : character(len=128) :: errmsg
998 : character(len=*), parameter :: subname = 'get_chunk_info_p: '
999 :
1000 1944000 : if (.not. phys_grid_initialized()) then
1001 0 : call endrun(subname//': physics grid not initialized')
1002 1944000 : else if ((index < 1) .or. (index > columns_on_task)) then
1003 0 : write(errmsg, '(a,2(a,i0))') subname, 'index (', index, &
1004 0 : ') out of range (1 to ', columns_on_task
1005 0 : write(iulog, *) errmsg
1006 0 : call endrun(errmsg)
1007 : else
1008 1944000 : lchnk = phys_columns(index)%local_phys_chunk
1009 1944000 : icol = phys_columns(index)%phys_chunk_index
1010 : end if
1011 :
1012 1944000 : end subroutine get_chunk_info_p
1013 :
1014 : !========================================================================
1015 :
1016 0 : subroutine get_grid_dims(hdim1_d_out, hdim2_d_out)
1017 : use cam_abortutils, only: endrun
1018 : ! retrieve dynamics field grid information
1019 : ! hdim1_d and hdim2_d are dimensions of rectangular horizontal grid
1020 : ! data structure, If 1D data structure, then hdim2_d == 1.
1021 : integer, intent(out) :: hdim1_d_out
1022 : integer, intent(out) :: hdim2_d_out
1023 :
1024 0 : if (.not. phys_grid_initialized()) then
1025 0 : call endrun('get_grid_dims: physics grid not initialized')
1026 : end if
1027 0 : hdim1_d_out = hdim1_d
1028 0 : hdim2_d_out = hdim2_d
1029 :
1030 0 : end subroutine get_grid_dims
1031 :
1032 : !========================================================================
1033 :
1034 : ! Note: This routine is a stub for future load-balancing
1035 0 : subroutine phys_decomp_to_dyn()
1036 : !-----------------------------------------------------------------------
1037 : !
1038 : ! phys_decomp_to_dyn: Transfer physics data to dynamics decomp
1039 : !
1040 : !-----------------------------------------------------------------------
1041 0 : end subroutine phys_decomp_to_dyn
1042 :
1043 : !========================================================================
1044 :
1045 : ! Note: This routine is a stub for future load-balancing
1046 0 : subroutine dyn_decomp_to_phys()
1047 : !-----------------------------------------------------------------------
1048 : !
1049 : ! dyn_decomp_to_phys: Transfer dynamics data to physics decomp
1050 : !
1051 : !-----------------------------------------------------------------------
1052 :
1053 0 : end subroutine dyn_decomp_to_phys
1054 :
1055 : !========================================================================
1056 :
1057 : subroutine dump_grid_map(grid_map)
1058 : use spmd_utils, only: iam, npes, mpicom
1059 : use cam_grid_support, only: iMap
1060 :
1061 : integer(iMap), pointer :: grid_map(:,:)
1062 :
1063 : integer :: num_cols
1064 : integer :: penum, icol
1065 : logical :: unstruct
1066 : integer :: file
1067 : integer :: ierr
1068 :
1069 : unstruct = SIZE(grid_map, 1) == 3
1070 : num_cols = SIZE(grid_map, 2)
1071 : if (iam == 0) then
1072 : open(newunit=file, file='physgrid_map.csv', status='replace')
1073 : if (unstruct) then
1074 : write(file, *) '"iam","col","block","map pos"'
1075 : else
1076 : write(file, *) '"iam","col","block","lon","lat"'
1077 : end if
1078 : close(unit=file)
1079 : end if
1080 : do penum = 0, npes - 1
1081 : if (iam == penum) then
1082 : open(newunit=file, file='physgrid_map.csv', status='old', &
1083 : action='readwrite', position='append')
1084 : do icol = 1, num_cols
1085 : if (unstruct) then
1086 : write(file, '(3(i0,","),i0)') iam, int(grid_map(1,icol)), &
1087 : int(grid_map(2,icol)), int(grid_map(3,icol))
1088 : else
1089 : write(file, '(4(i0,","),i0)') iam, int(grid_map(1,icol)), &
1090 : int(grid_map(2,icol)), int(grid_map(3,icol)), &
1091 : int(grid_map(4,icol))
1092 : end if
1093 : end do
1094 : close(unit=file)
1095 : end if
1096 : call MPI_barrier(mpicom, ierr)
1097 : end do
1098 : end subroutine dump_grid_map
1099 :
1100 : !=============================================================================
1101 : !==
1102 : !!!!!! DUMMY INTERFACEs TO TEST WEAK SCALING INFRASTRUCTURE, SHOULD GO AWAY
1103 : !==
1104 : !=============================================================================
1105 :
1106 0 : subroutine scatter_field_to_chunk(fdim,mdim,ldim, &
1107 : hdim1d,globalfield,localchunks)
1108 : use cam_abortutils, only: endrun
1109 : use ppgrid, only: pcols
1110 : !-----------------------------------------------------------------------
1111 : !
1112 : ! Purpose: DUMMY FOR WEAK SCALING TESTS
1113 : !
1114 : !------------------------------Arguments--------------------------------
1115 : integer, intent(in) :: fdim ! declared length of first dimension
1116 : integer, intent(in) :: mdim ! declared length of middle dimension
1117 : integer, intent(in) :: ldim ! declared length of last dimension
1118 : integer, intent(in) :: hdim1d ! declared first horizontal index
1119 : real(r8), intent(in) :: globalfield(fdim,hdim1d,mdim,hdim2_d,ldim)
1120 : real(r8), intent(out):: localchunks(fdim,pcols,mdim, &
1121 : begchunk:endchunk,ldim)
1122 :
1123 0 : call endrun('scatter_field_to_chunk: NOT SUPPORTED WITH WEAK SCALING')
1124 0 : end subroutine scatter_field_to_chunk
1125 :
1126 : !========================================================================
1127 :
1128 0 : subroutine get_lat_all_p_int(lcid, latdim, lats)
1129 : use cam_abortutils, only: endrun
1130 : !-----------------------------------------------------------------------
1131 : !
1132 : ! get_lat_all_p: Return all latitudes (in degrees) for chunk, <lcid>
1133 : !
1134 : !-----------------------------------------------------------------------
1135 : ! Dummy Arguments
1136 : integer, intent(in) :: lcid ! local chunk id
1137 : integer, intent(in) :: latdim ! declared size of output array
1138 : integer, intent(out) :: lats(latdim) ! array of latitudes
1139 :
1140 0 : call endrun('get_lat_all_p: deprecated interface')
1141 :
1142 0 : end subroutine get_lat_all_p_int
1143 :
1144 : !========================================================================
1145 :
1146 0 : subroutine get_lon_all_p_int(lcid, londim, lons)
1147 : use cam_abortutils, only: endrun
1148 : !-----------------------------------------------------------------------
1149 : !
1150 : ! get_lon_all_p:: Return all longitudes (in degrees) for chunk, <lcid>
1151 : !
1152 : !-----------------------------------------------------------------------
1153 : ! Dummy Arguments
1154 : integer, intent(in) :: lcid ! local chunk id
1155 : integer, intent(in) :: londim ! declared size of output array
1156 : integer, intent(out) :: lons(londim) ! array of longitudes
1157 :
1158 0 : call endrun('get_lon_all_p: deprecated interface')
1159 :
1160 0 : end subroutine get_lon_all_p_int
1161 :
1162 : !========================================================================
1163 :
1164 0 : end module phys_grid
|