MPAS-JEDI
mpas_geom_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2017 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
8 use atlas_module, only: atlas_functionspace, atlas_fieldset, atlas_field, atlas_real
9 use fckit_configuration_module, only: fckit_configuration, fckit_yamlconfiguration
10 use fckit_pathname_module, only: fckit_pathname
11 use fckit_log_module, only: fckit_log
12 use fckit_mpi_module, only: fckit_mpi_comm
13 use iso_c_binding
14 
15 !oops
16 use oops_variables_mod, only: oops_variables
17 
18 !ufo
19 use ufo_vars_mod, only: maxvarlen, ufo_vars_getindex
20 
21 !MPAS-Model
22 use mpas_derived_types
23 use mpas_kind_types
24 use mpas_constants
25 use kinds, only : kind_real
26 use mpas_dmpar, only: mpas_dmpar_sum_int
27 use mpas_subdriver
28 use atm_core
29 use mpas_pool_routines
30 
31 !mpas_jedi
33 
34 implicit none
35 private
36 public :: mpas_geom, &
40 
41 public :: mpas_geom_registry
42 
43 
44 ! ------------------------------------------------------------------------------
45 
47  character(len=MAXVARLEN) :: name
48  character(len=MAXVARLEN) :: template
49  character(len=MAXVARLEN) :: identity
50 end type templated_field
51 
52 !> Fortran derived type to hold geometry definition
53 type :: mpas_geom
54  integer :: ncellsglobal !Global count
55  integer :: nedgesglobal !Global count
56  integer :: nverticesglobal !Global count
57  integer :: ncells !Memory count (Local + Halo)
58  integer :: nedges !Memory count (Local + Halo)
59  integer :: nvertices !Memory count (Local + Halo)
60  integer :: ncellssolve !Local count
61  integer :: nedgessolve !Local count
62  integer :: nverticessolve !Local count
63  integer :: nvertlevels
64  integer :: nvertlevelsp1
65  integer :: nsoillevels
66  integer :: vertexdegree
67  integer :: maxedges
68  logical :: deallocate_nonda_fields
69  logical :: use_bump_interpolation
70  character(len=StrKIND) :: bump_vunit
71  real(kind=kind_real), dimension(:), allocatable :: latcell, loncell
72  real(kind=kind_real), dimension(:), allocatable :: areacell
73  real(kind=kind_real), dimension(:), allocatable :: latedge, lonedge
74  real(kind=kind_real), dimension(:,:), allocatable :: edgenormalvectors
75  real(kind=kind_real), dimension(:,:), allocatable :: zgrid
76  integer, allocatable :: nedgesoncell(:)
77  integer, allocatable :: cellsoncell(:,:)
78  integer, allocatable :: edgesoncell(:,:)
79  integer, allocatable :: cellsonvertex(:,:)
80  integer, allocatable :: cellsonedge(:,:)
81  integer, allocatable :: verticesonedge(:,:)
82  real(kind=kind_real), DIMENSION(:), ALLOCATABLE :: dcedge, dvedge
83  real(kind=kind_real), DIMENSION(:), ALLOCATABLE :: areatriangle, angleedge
84  real(kind=kind_real), DIMENSION(:,:), ALLOCATABLE :: kiteareasonvertex, edgesoncell_sign
85 
86  type (domain_type), pointer :: domain => null()
87  type (core_type), pointer :: corelist => null()
88 
89  type(fckit_mpi_comm) :: f_comm
90 
91  type(atlas_functionspace) :: afunctionspace
92 
93  type(templated_field), allocatable :: templated_fields(:)
94 
95  contains
96 
97  procedure, public :: is_templated => field_is_templated
98  procedure, public :: template => template_fieldname
99  procedure, public :: has_identity => field_has_identity
100  procedure, public :: identity => identity_fieldname
101  generic, public :: nlevels => variables_nlevels, var_nlevels
102  procedure :: variables_nlevels
103  procedure :: var_nlevels
104 
105 end type mpas_geom
106 
107 type :: idcounter
108  integer :: id, counter
109 end type idcounter
110 
111 type(idcounter), allocatable :: geom_count(:)
112 
113 character(len=1024) :: message
114 
115 character(len=MAXVARLEN), parameter :: &
117  [character(len=maxvarlen) :: 'nVertLevels', 'nVertLevelsP1', &
118  'nVertLevelsP2', 'nSoilLevels', 'nOznLevels', 'nAerLevels']
119 
120 #define LISTED_TYPE mpas_geom
121 
122 !> Linked list interface - defines registry_t type
123 #include <oops/util/linkedList_i.f>
124 
125 !> Global registry
126 type(registry_t) :: mpas_geom_registry
127 
128 ! ------------------------------------------------------------------------------
129 contains
130 ! ------------------------------------------------------------------------------
131 !> Linked list implementation
132 #include <oops/util/linkedList_c.f>
133 
134 ! ------------------------------------------------------------------------------
135 subroutine geo_setup(self, f_conf, f_comm)
136 
137  implicit none
138 
139  type(mpas_geom), intent(inout) :: self
140  type(fckit_configuration), intent(in) :: f_conf
141  type(fckit_mpi_comm), intent(in) :: f_comm
142 
143  character(len=StrKIND) :: string1
144  type (mpas_pool_type), pointer :: meshpool, fg
145  type (block_type), pointer :: block_ptr
146  !character(len=120) :: fn
147  character(len=512) :: nml_file, streams_file, fields_file
148  character(len=:), allocatable :: str
149  type(fckit_configuration) :: template_conf
150  type(fckit_configuration), allocatable :: fields_conf(:)
151 
152  real (kind=kind_real), pointer :: r1d_ptr(:), r2d_ptr(:,:)
153  integer, pointer :: i0d_ptr, i1d_ptr(:), i2d_ptr(:,:)
154 
155  type(idcounter), allocatable :: prev_count(:)
156  integer :: ii, nprev
157 
158  logical :: deallocate_fields
159  logical :: bump_interp
160 
161  call fckit_log%info('==> create geom')
162 
163  ! MPI communicator
164  self%f_comm = f_comm
165 
166  ! "nml_file" and "streams_file" are mandatory for mpas_init.
167  call f_conf%get_or_die("nml_file",str)
168  nml_file = str
169  call f_conf%get_or_die("streams_file",str)
170  streams_file = str
171 
172  ! Domain decomposition and templates for state/increment variables
173  call mpas_init( self % corelist, self % domain, mpi_comm = self%f_comm%communicator(), &
174  & namelistfileparam = trim(nml_file), streamsfileparam = trim(streams_file))
175 
176  ! The interpolation method specified here will be used to interpolate data between
177  ! different geometries, but not in GetValues. GetValues has its own configuration
178  ! for specifying its interpolation method.
179  if (f_conf%get("interpolation type", str)) then
180  select case (str)
181  case ('bump')
182  self%use_bump_interpolation = .true.
183  case ('unstructured')
184  self%use_bump_interpolation = .false.
185  case default
186  write(message,*) '--> geo_setup: interpolation type: ',str,' not implemented'
187  call abor1_ftn(message)
188  end select
189  else
190  self%use_bump_interpolation = .true. ! BUMP is default interpolation
191  end if
192 
193  !Deallocate not-used fields for memory reduction
194  if (f_conf%has("deallocate non-da fields")) then
195  call f_conf%get_or_die("deallocate non-da fields",deallocate_fields)
196  self % deallocate_nonda_fields = deallocate_fields
197  else
198  self % deallocate_nonda_fields = .false.
199  end if
200  if (self % deallocate_nonda_fields) call geo_deallocate_nonda_fields (self % domain)
201 
202  ! Set up the vertical coordinate for bump
203  if (f_conf%has("bump vunit")) then
204  call f_conf%get_or_die("bump vunit",str)
205  self % bump_vunit = str
206  else
207  self % bump_vunit = 'modellevel'
208  end if
209 
210  if (allocated(geom_count)) then
211  nprev = size(geom_count)
212  allocate(prev_count(nprev))
213  do ii = 1, nprev
214  prev_count(ii) = geom_count(ii)
215  if (prev_count(ii)%id == self%domain%domainID) then
216  call abor1_ftn("domainID already used")
217  end if
218  end do
219  deallocate(geom_count)
220  else
221  nprev = 0
222  end if
223  allocate(geom_count(nprev+1))
224  do ii = 1, nprev
225  geom_count(ii) = prev_count(ii)
226  end do
227  geom_count(nprev+1)%id = self%domain%domainID
228  geom_count(nprev+1)%counter = 1
229 
230  ! first read array of templated field configurations
231  if (f_conf%get('template fields file',str)) then
232  fields_file = str
233  else
234  fields_file = 'geovars.yaml'
235  end if
236  template_conf = fckit_yamlconfiguration(fckit_pathname(fields_file))
237  call template_conf%get_or_die('fields',fields_conf)
238 
239  ! create templated field descriptors
240  allocate(self%templated_fields(size(fields_conf)))
241  do ii = 1, size(fields_conf)
242  call fields_conf(ii)%get_or_die('field name',str)
243  self%templated_fields(ii)%name = trim(str)
244  call fields_conf(ii)%get_or_die('mpas template field',str)
245  self%templated_fields(ii)%template = trim(str)
246  if (fields_conf(ii)%get('mpas identity field',str)) then
247  self%templated_fields(ii)%identity = trim(str)
248  else
249  self%templated_fields(ii)%identity = 'none'
250  end if
251  end do
252  deallocate(fields_conf)
253 
254 ! if (associated(self % domain)) then
255 ! write(message,*) 'inside geom: geom % domain associated for domainID = ', self % domain % domainID
256 ! call fckit_log%debug(message)
257 ! end if
258 ! if (associated(self % corelist)) then
259 ! call fckit_log%debug('inside geom: geom % corelist associated')
260 ! else
261 ! call fckit_log%debug('inside geom: geom % corelist not associated')
262 ! end if
263 
264  ! These pool accesses refer to memory (local+halo) for a single MPAS block (standard)
265  block_ptr => self % domain % blocklist
266 
267  call mpas_pool_get_subpool ( block_ptr % structs, 'mesh', meshpool )
268 
269  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nCells', i0d_ptr )
270  self % nCells = i0d_ptr
271  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nCellsSolve', i0d_ptr )
272  self % nCellsSolve = i0d_ptr
273  call mpas_dmpar_sum_int ( self % domain % dminfo, &
274  self % nCellsSolve, self % nCellsGlobal )
275 
276  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nEdges', i0d_ptr )
277  self % nEdges = i0d_ptr
278  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nEdgesSolve', i0d_ptr )
279  self % nEdgesSolve = i0d_ptr
280  call mpas_dmpar_sum_int ( self % domain % dminfo, &
281  self % nEdgesSolve, self % nEdgesGlobal )
282 
283  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nVertices', i0d_ptr )
284  self % nVertices = i0d_ptr
285  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nVerticesSolve', i0d_ptr )
286  self % nVerticesSolve = i0d_ptr
287  call mpas_dmpar_sum_int ( self % domain % dminfo, &
288  self % nVerticesSolve, self % nVerticesGlobal )
289 
290  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nVertLevels', i0d_ptr )
291  self % nVertLevels = i0d_ptr
292  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nVertLevelsP1', i0d_ptr )
293  self % nVertLevelsP1 = i0d_ptr
294  call mpas_pool_get_dimension ( block_ptr % dimensions, 'nSoilLevels', i0d_ptr )
295  self % nSoilLevels = i0d_ptr
296  call mpas_pool_get_dimension ( block_ptr % dimensions, 'vertexDegree', i0d_ptr )
297  self % vertexDegree = i0d_ptr
298  call mpas_pool_get_dimension ( block_ptr % dimensions, 'maxEdges', i0d_ptr )
299  self % maxEdges = i0d_ptr
300 
301  allocate ( self % latCell ( self % nCells ) )
302  allocate ( self % lonCell ( self % nCells ) )
303  allocate ( self % latEdge ( self % nEdges ) )
304  allocate ( self % lonEdge ( self % nEdges ) )
305  allocate ( self % areaCell ( self % nCells ) )
306  allocate ( self % edgeNormalVectors (3, self % nEdges ) )
307  allocate ( self % zgrid ( self % nVertLevelsP1, self % nCells ) )
308  allocate ( self % nEdgesOnCell ( self % nCells ) )
309  allocate ( self % edgesOnCell ( self % maxEdges, self % nCells ) )
310  allocate ( self % cellsOnCell ( self % maxEdges, self % nCells ) )
311  allocate ( self % cellsOnVertex ( self % vertexDegree, self % nVertices ) )
312  allocate ( self % cellsOnEdge ( 2, self % nEdges ) )
313  allocate ( self % verticesOnEdge ( 2, self % nEdges ) )
314  allocate ( self % dcEdge ( self % nEdges ) )
315  allocate ( self % dvEdge ( self % nEdges ) )
316  allocate ( self % kiteAreasOnVertex ( self % vertexDegree, self % nVertices ) )
317  allocate ( self % edgesOnCell_sign ( self % maxEdges, self % nCells ) )
318  allocate ( self % areaTriangle ( self % nVertices ) )
319  allocate ( self % angleEdge ( self % nEdges ) )
320 
321  call mpas_pool_get_array ( meshpool, 'latCell', r1d_ptr )
322  self % latCell = r1d_ptr(1:self % nCells)
323  where (self % latCell > mpas_jedi_piio2_kr)
324  self % latCell = mpas_jedi_piio2_kr
325  end where
326  where (self % latCell < - mpas_jedi_piio2_kr)
327  self % latCell = - mpas_jedi_piio2_kr
328  end where
329 
330  call mpas_pool_get_array ( meshpool, 'lonCell', r1d_ptr )
331  self % lonCell = r1d_ptr(1:self % nCells)
332  call mpas_pool_get_array ( meshpool, 'areaCell', r1d_ptr )
333  self % areaCell = r1d_ptr(1:self % nCells)
334  call mpas_pool_get_array ( meshpool, 'latEdge', r1d_ptr )
335  self % latEdge = r1d_ptr(1:self % nEdges)
336  call mpas_pool_get_array ( meshpool, 'lonEdge', r1d_ptr )
337  self % lonEdge = r1d_ptr(1:self % nEdges)
338  call mpas_pool_get_array ( meshpool, 'edgeNormalVectors', r2d_ptr )
339  self % edgeNormalVectors = r2d_ptr( 1:3, 1:self % nEdges )
340  call mpas_pool_get_array ( meshpool, 'nEdgesOnCell', i1d_ptr )
341  self % nEdgesOnCell = i1d_ptr(1:self % nCells)
342  call mpas_pool_get_array ( meshpool, 'edgesOnCell', i2d_ptr )
343  self % edgesOnCell = i2d_ptr( 1:self % maxEdges, 1:self % nCells )
344  call mpas_pool_get_array ( meshpool, 'cellsOnCell', i2d_ptr )
345  self % cellsOnCell = i2d_ptr( 1:self % maxEdges, 1:self % nCells )
346  call mpas_pool_get_array ( meshpool, 'cellsOnVertex', i2d_ptr )
347  self % cellsOnVertex = i2d_ptr( 1:self % vertexDegree, 1:self % nVertices )
348  call mpas_pool_get_array ( meshpool, 'cellsOnEdge', i2d_ptr )
349  self % cellsOnEdge = i2d_ptr( 1:2, 1:self % nEdges )
350  call mpas_pool_get_array ( meshpool, 'verticesOnEdge', i2d_ptr )
351  self % verticesOnEdge = i2d_ptr( 1:2, 1:self % nEdges )
352  call mpas_pool_get_array ( meshpool, 'dcEdge', r1d_ptr )
353  self % dcEdge = r1d_ptr(1:self % nEdges)
354  call mpas_pool_get_array ( meshpool, 'dvEdge', r1d_ptr )
355  self % dvEdge = r1d_ptr(1:self % nEdges)
356  call mpas_pool_get_array ( meshpool, 'kiteAreasOnVertex', r2d_ptr )
357  self % kiteAreasOnVertex = r2d_ptr( 1:self % vertexDegree, 1:self % nVertices )
358  call mpas_pool_get_array ( meshpool, 'edgesOnCell_sign', r2d_ptr )
359  self % edgesOnCell_sign = r2d_ptr( 1:self % maxEdges, 1:self % nCells )
360  call mpas_pool_get_array ( meshpool, 'areaTriangle', r1d_ptr )
361  self % areaTriangle = r1d_ptr(1:self % nVertices)
362  call mpas_pool_get_array ( meshpool, 'angleEdge', r1d_ptr )
363  self % angleEdge = r1d_ptr(1:self % nEdges)
364 
365  call mpas_pool_get_array ( meshpool, 'zgrid', r2d_ptr )
366  self % zgrid = r2d_ptr( 1:self % nVertLevelsP1, 1:self % nCells )
367 
368  call fckit_log%debug('End of geo_setup')
369  if (allocated(prev_count)) deallocate(prev_count)
370  if (allocated(str)) deallocate(str)
371 
372 end subroutine geo_setup
373 
374 ! --------------------------------------------------------------------------------------------------
375 
376 subroutine geo_set_atlas_lonlat(self, afieldset)
377 
378  implicit none
379 
380  type(mpas_geom), intent(inout) :: self
381  type(atlas_fieldset), intent(inout) :: afieldset
382 
383  real(kind_real), pointer :: real_ptr(:,:)
384  type(atlas_field) :: afield
385 
386  ! Create lon/lat field
387  afield = atlas_field(name="lonlat", kind=atlas_real(kind_real), shape=(/2,self%nCellsSolve/))
388  call afield%data(real_ptr)
389  real_ptr(1,:) = self%lonCell(1:self%nCellsSolve) * mpas_jedi_rad2deg_kr
390  real_ptr(2,:) = self%latCell(1:self%nCellsSolve) * mpas_jedi_rad2deg_kr
391  call afieldset%add(afield)
392 
393 end subroutine geo_set_atlas_lonlat
394 
395 ! --------------------------------------------------------------------------------------------------
396 
397 subroutine geo_fill_atlas_fieldset(self, afieldset)
398 
399  implicit none
400 
401  type(mpas_geom), intent(inout) :: self
402  type(atlas_fieldset), intent(inout) :: afieldset
403 
404  integer :: i, jz
405  real(kind=kind_real), pointer :: real_ptr_1(:), real_ptr_2(:,:)
406  type(atlas_field) :: afield
407 
408  ! Add area
409  afield = self%afunctionspace%create_field(name='area', kind=atlas_real(kind_real), levels=0)
410  call afield%data(real_ptr_1)
411  real_ptr_1 = self%areaCell(1:self%nCellsSolve)
412  call afieldset%add(afield)
413  call afield%final()
414 
415  ! Add vertical unit
416  afield = self%afunctionspace%create_field(name='vunit', kind=atlas_real(kind_real), levels=self%nVertLevels)
417  call afield%data(real_ptr_2)
418  do jz=1,self%nVertLevels
419  if (trim(self % bump_vunit) .eq. 'modellevel') then
420  real_ptr_2(jz,:) = real(jz, kind_real)
421  else if (trim(self % bump_vunit) .eq. 'height') then
422  real_ptr_2(jz,1:self%nCellsSolve) = mpas_jedi_half_kr * &
423  ( self%zgrid(jz,1:self%nCellsSolve) + self%zgrid(jz+1,1:self%nCellsSolve) )
424  end if
425  !--Similarly for ln of pressure_base. I don't know if we can access to "total pressure" here.
426  !real (kind=kind_real), dimension(:,:), pointer :: pressure_base
427  !call mpas_pool_get_array(self%domain%blocklist%allFields,'pressure_base', pressure_base)
428  !real_ptr_2(jz,1:self%nCellsSolve) = log( pressure_base(jz,1:self%nCellsSolve) )
429  end do
430  call afieldset%add(afield)
431  call afield%final()
432 
433 end subroutine geo_fill_atlas_fieldset
434 
435 ! ------------------------------------------------------------------------------
436 subroutine geo_deallocate_nonda_fields(domain)
437 
438  implicit none
439  type (domain_type), pointer, intent(inout) :: domain
440  type (mpas_pool_type), pointer :: pool_a, pool_b
441  type (mpas_pool_data_type), pointer :: mem
442  type (mpas_pool_iterator_type) :: poolItr_b
443 
444  type (field0DReal), pointer :: field0d
445  type (field1DReal), pointer :: field1d
446  type (field2DReal), pointer :: field2d
447 
448  integer, parameter :: num_da_fields = 36
449  integer :: i
450  character (len=22), allocatable :: poolname_a(:), poolname_b(:)
451  character (len=22), allocatable :: da_fieldnames(:)
452 
453  allocate(poolname_a(3))
454  poolname_a(1)='tend'
455  poolname_a(2)='tend_physics'
456  poolname_a(3)='atm_input'
457 
458  allocate(poolname_b(3))
459  poolname_b(1)='diag_physics'
460  poolname_b(2)='diag'
461  poolname_b(3)='sfc_input'
462 
463  !TODO: remove the following list of fieldnames and read from a stream_list.atmosphere file
464  allocate(da_fieldnames(num_da_fields))
465  da_fieldnames(1)='re_cloud'
466  da_fieldnames(2)='re_ice'
467  da_fieldnames(3)='re_snow'
468  da_fieldnames(4)='cldfrac'
469  da_fieldnames(5)='lai'
470  da_fieldnames(6)='u10'
471  da_fieldnames(7)='v10'
472  da_fieldnames(8)='q2'
473  da_fieldnames(9)='t2m'
474  da_fieldnames(10)='pressure_p'
475  da_fieldnames(11)='pressure'
476  da_fieldnames(12)='pressure_base'
477  da_fieldnames(13)='surface_pressure'
478  da_fieldnames(14)='rho'
479  da_fieldnames(15)='theta'
480  da_fieldnames(16)='temperature'
481  da_fieldnames(17)='relhum'
482  da_fieldnames(18)='spechum'
483  da_fieldnames(19)='u'
484  da_fieldnames(20)='uReconstructZonal'
485  da_fieldnames(21)='uReconstructMeridional'
486  da_fieldnames(22)='stream_function'
487  da_fieldnames(23)='velocity_potential'
488  da_fieldnames(24)='landmask'
489  da_fieldnames(25)='xice'
490  da_fieldnames(26)='snowc'
491  da_fieldnames(27)='skintemp'
492  da_fieldnames(28)='ivgtyp'
493  da_fieldnames(29)='isltyp'
494  da_fieldnames(30)='snowh'
495  da_fieldnames(31)='vegfra'
496  da_fieldnames(32)='lai'
497  da_fieldnames(33)='smois'
498  da_fieldnames(34)='tslb'
499  da_fieldnames(35)='vorticity'
500  da_fieldnames(36)='w'
501 
502  do i=1,size(poolname_a)
503  mem => pool_get_member(domain % blocklist % structs, poolname_a(i), mpas_pool_subpool)
504  if (associated(mem)) then
505  call mpas_pool_get_subpool(domain % blocklist % structs, poolname_a(i), pool_a)
506  call mpas_pool_destroy_pool(pool_a)
507  call mpas_pool_remove_subpool(domain % blocklist % structs, poolname_a(i))
508  end if
509  end do
510 
511  do i=1,size(poolname_b)
512 
513  mem => pool_get_member(domain % blocklist % structs, poolname_b(i), mpas_pool_subpool)
514  if (associated(mem)) then
515  call mpas_pool_get_subpool(domain % blocklist % structs, poolname_b(i), pool_b)
516  call mpas_pool_begin_iteration(pool_b)
517 
518  do while ( mpas_pool_get_next_member(pool_b, poolitr_b) )
519 
520  if (poolitr_b % memberType == mpas_pool_field .and. poolitr_b % dataType == mpas_pool_real .and. &
521  ufo_vars_getindex(da_fieldnames, poolitr_b % memberName) < 1) then
522 
523  if (poolitr_b % nDims == 0) then
524  call mpas_pool_get_field(pool_b,trim(poolitr_b % memberName),field0d)
525  if (associated(field0d)) then
526  call mpas_deallocate_field(field0d)
527  call mpas_pool_remove_field(pool_b, trim(poolitr_b % memberName))
528  end if
529  nullify(field0d)
530  else if (poolitr_b % nDims == 1) then
531  call mpas_pool_get_field(pool_b,trim(poolitr_b % memberName),field1d)
532  if (associated(field1d)) then
533  call mpas_deallocate_field(field1d)
534  call mpas_pool_remove_field(pool_b, trim(poolitr_b % memberName))
535  end if
536  nullify(field1d)
537  else if (poolitr_b % nDims == 2) then
538  call mpas_pool_get_field(pool_b,trim(poolitr_b % memberName),field2d)
539  if (associated(field2d)) then
540  call mpas_deallocate_field(field2d)
541  call mpas_pool_remove_field(pool_b, trim(poolitr_b % memberName))
542  end if
543  nullify(field2d)
544  endif
545 
546  end if
547 
548  end do
549 
550  end if
551 
552  end do
553 
554  deallocate(poolname_a)
555  deallocate(poolname_b)
556  deallocate(da_fieldnames)
557 
558 end subroutine geo_deallocate_nonda_fields
559 
560 ! ------------------------------------------------------------------------------
561 
562 subroutine geo_clone(self, other)
563 
564  implicit none
565 
566  type(mpas_geom), intent(inout) :: self
567  type(mpas_geom), intent(in) :: other
568 
569  integer :: ii
570 
571  ! Clone communicator
572  self%f_comm = other%f_comm
573 
574  call fckit_log%debug('====> copy of geom array')
575 
576  self % nCellsGlobal = other % nCellsGlobal
577  self % nCells = other % nCells
578  self % nCellsSolve = other % nCellsSolve
579 
580  self % nEdgesGlobal = other % nEdgesGlobal
581  self % nEdges = other % nEdges
582  self % nEdgesSolve = other % nEdgesSolve
583 
584  self % nVerticesGlobal = other % nVerticesGlobal
585  self % nVertices = other % nVertices
586  self % nVerticesSolve = other % nVerticesSolve
587 
588  self % nVertLevels = other % nVertLevels
589  self % nVertLevelsP1 = other % nVertLevelsP1
590  self % nSoilLevels = other % nSoilLevels
591  self % vertexDegree = other % vertexDegree
592  self % maxEdges = other % maxEdges
593 
594  if (.not.allocated(self % latCell)) allocate(self % latCell(other % nCells))
595  if (.not.allocated(self % lonCell)) allocate(self % lonCell(other % nCells))
596  if (.not.allocated(self % latEdge)) allocate(self % latEdge(other % nEdges))
597  if (.not.allocated(self % lonEdge)) allocate(self % lonEdge(other % nEdges))
598  if (.not.allocated(self % areaCell)) allocate(self % areaCell(other % nCells))
599  if (.not.allocated(self % edgeNormalVectors)) allocate(self % edgeNormalVectors(3, other % nEdges))
600  if (.not.allocated(self % zgrid)) allocate(self % zgrid(other % nVertLevelsP1, other % nCells))
601  if (.not.allocated(self % nEdgesOnCell)) allocate(self % nEdgesOnCell(other % nCells))
602  if (.not.allocated(self % edgesOnCell)) allocate(self % edgesOnCell(other % maxEdges, other % nCells))
603  if (.not.allocated(self % cellsOnCell)) allocate(self % cellsOnCell (other % maxEdges, other % nCells))
604  if (.not.allocated(self % cellsOnVertex)) allocate(self % cellsOnVertex(self % vertexDegree, self % nVertices))
605  if (.not.allocated(self % cellsOnEdge)) allocate(self % cellsOnEdge (2, self % nEdges))
606  if (.not.allocated(self % verticesOnEdge)) allocate(self % verticesOnEdge (2, self % nEdges))
607  if (.not.allocated(self % dcEdge)) allocate(self % dcEdge (self % nEdges))
608  if (.not.allocated(self % dvEdge)) allocate(self % dvEdge (self % nEdges))
609  if (.not.allocated(self % kiteAreasOnVertex)) allocate(self % kiteAreasOnVertex(self % vertexDegree, self % nVertices))
610  if (.not.allocated(self % edgesOnCell_sign)) allocate(self % edgesOnCell_sign(self % maxEdges, self % nCells))
611  if (.not.allocated(self % areaTriangle)) allocate(self % areaTriangle(self % nVertices))
612  if (.not.allocated(self % angleEdge)) allocate(self % angleEdge(self % nEdges))
613 
614  self % use_bump_interpolation = other % use_bump_interpolation
615  self % templated_fields = other % templated_fields
616  self % latCell = other % latCell
617  self % lonCell = other % lonCell
618  self % areaCell = other % areaCell
619  self % latEdge = other % latEdge
620  self % lonEdge = other % lonEdge
621  self % edgeNormalVectors = other % edgeNormalVectors
622  self % zgrid = other % zgrid
623  self % nEdgesOnCell = other % nEdgesOnCell
624  self % edgesOnCell = other % edgesOnCell
625  self % cellsOnCell = other % cellsOnCell
626  self % cellsOnVertex = other % cellsOnVertex
627  self % cellsOnEdge = other % cellsOnEdge
628  self % verticesOnEdge = other % verticesOnEdge
629  self % dcEdge = other % dcEdge
630  self % dvEdge = other % dvEdge
631  self % kiteAreasOnVertex = other % kiteAreasOnVertex
632  self % edgesOnCell_sign = other % edgesOnCell_sign
633  self % areaTriangle = other % areaTriangle
634  self % angleEdge = other % angleEdge
635 
636  call fckit_log%debug('====> copy of geom corelist and domain')
637 
638  self % corelist => other % corelist
639  self % domain => other % domain
640  write(message,*) 'inside geo_clone: other % domain % domainID = ', other % domain % domainID
641  call fckit_log%debug(message)
642 
643  do ii = 1, size(geom_count)
644  if (geom_count(ii)%id == self%domain%domainID) then
645  geom_count(ii)%counter = geom_count(ii)%counter + 1
646  end if
647  end do
648 
649  call fckit_log%debug('====> copy of geom done')
650 
651 end subroutine geo_clone
652 
653 ! ------------------------------------------------------------------------------
654 
655 subroutine geo_delete(self)
656 
657  implicit none
658 
659  type(mpas_geom), intent(inout) :: self
660  integer :: ii
661 
662  if (allocated(self%templated_fields)) deallocate(self%templated_fields)
663  if (allocated(self%latCell)) deallocate(self%latCell)
664  if (allocated(self%lonCell)) deallocate(self%lonCell)
665  if (allocated(self%latEdge)) deallocate(self%latEdge)
666  if (allocated(self%lonEdge)) deallocate(self%lonEdge)
667  if (allocated(self%areaCell)) deallocate(self%areaCell)
668  if (allocated(self%edgeNormalVectors)) deallocate(self%edgeNormalVectors)
669  if (allocated(self%zgrid)) deallocate(self%zgrid)
670  if (allocated(self%nEdgesOnCell)) deallocate(self%nEdgesOnCell)
671  if (allocated(self%edgesOnCell)) deallocate(self%edgesOnCell)
672  if (allocated(self%cellsOnCell)) deallocate(self%cellsOnCell)
673  if (allocated(self%cellsOnVertex)) deallocate(self%cellsOnVertex)
674  if (allocated(self%cellsOnEdge)) deallocate(self%cellsOnEdge)
675  if (allocated(self%verticesOnEdge)) deallocate(self%verticesOnEdge)
676  if (allocated(self%dcEdge)) deallocate(self%dcEdge)
677  if (allocated(self%dvEdge)) deallocate(self%dvEdge)
678  if (allocated(self%kiteAreasOnVertex)) deallocate(self%kiteAreasOnVertex)
679  if (allocated(self%edgesOnCell_sign)) deallocate(self%edgesOnCell_sign)
680  if (allocated(self%areaTriangle)) deallocate(self%areaTriangle)
681  if (allocated(self%angleEdge)) deallocate(self%angleEdge)
682 
683  do ii = 1, size(geom_count)
684  if (geom_count(ii)%id == self%domain%domainID) then
685  geom_count(ii)%counter = geom_count(ii)%counter - 1
686  if ((associated(self % corelist)).and.(associated(self % domain))) then
687  if (geom_count(ii)%counter == 0) then
688  ! Completely destroy corelist and domain
689  ! + can only do this if they are not used by another copy of self
690  ! + otherwise the persistance of the underlying domain
691  ! and corelist is a known memory leak
692  write(message,'(A,I3)') &
693  '==> destruct MPAS corelist and domain:', self%domain%domainID
694  call fckit_log%info(message)
695  call mpas_timer_set_context( self % domain )
696  call mpas_finalize(self % corelist, self % domain)
697  else
698  nullify(self % corelist)
699  nullify(self % domain)
700  end if
701  exit
702  end if
703  end if
704  end do
705 
706 end subroutine geo_delete
707 
708 ! ------------------------------------------------------------------------------
709 
710 subroutine geo_is_equal(is_equal, self, other)
711 
712  implicit none
713 
714  logical*1, intent(inout) :: is_equal
715  type(mpas_geom), intent(in) :: self
716  type(mpas_geom), intent(in) :: other
717 
718  ! Only compares two of many attributes of the two geometries, but
719  ! sufficient for current needs.
720  is_equal = (self%nCellsGlobal .eq. other%nCellsGlobal) .and. &
721  (self%nVertLevels .eq. other%nVertLevels)
722 
723 end subroutine geo_is_equal
724 
725 ! ------------------------------------------------------------------------------
726 
727 subroutine geo_info(self, nCellsGlobal, nCells, nCellsSolve, &
728  nEdgesGlobal, nEdges, nEdgesSolve, &
729  nVertLevels, nVertLevelsP1)
730 
731  implicit none
732 
733  type(mpas_geom), intent(in) :: self
734  integer, intent(inout) :: ncellsglobal, ncells, ncellssolve
735  integer, intent(inout) :: nedgesglobal, nedges, nedgessolve
736  integer, intent(inout) :: nvertlevels
737  integer, intent(inout) :: nvertlevelsp1
738 
739  ncellsglobal = self%nCellsGlobal
740  ncells = self%nCells
741  ncellssolve = self%nCellsSolve
742 
743  nedgesglobal = self%nEdgesGlobal
744  nedges = self%nEdges
745  nedgessolve = self%nEdgesSolve
746 
747  nvertlevels = self%nVertLevels
748  nvertlevelsp1 = self%nVertLevelsP1
749 
750 end subroutine geo_info
751 
752 ! ------------------------------------------------------------------------------
753 
754 function field_is_templated(self, fieldname) result(is_templated)
755  class(mpas_geom), intent(in) :: self
756  character(len=*), intent(in) :: fieldname
757  integer :: ii
758  logical :: is_templated
759  is_templated = .false.
760  if (allocated(self%templated_fields)) then
761  do ii = 1, size(self%templated_fields)
762  if (trim(self%templated_fields(ii)%name) == trim(fieldname)) then
763  is_templated = .true.
764  return
765  end if
766  end do
767  end if
768 end function field_is_templated
769 
770 ! ------------------------------------------------------------------------------
771 
772 function template_fieldname(self, fieldname) result(template)
773  class(mpas_geom), intent(in) :: self
774  character(len=*), intent(in) :: fieldname
775  integer :: ii
776  character(len=MAXVARLEN) :: template
777  do ii = 1, size(self%templated_fields)
778  if (trim(self%templated_fields(ii)%name) == trim(fieldname)) then
779  template = self%templated_fields(ii)%template
780  return
781  end if
782  end do
783  write(message,*) '--> mpas_geom % template_fieldname: fieldname is not templated, ',fieldname
784  call abor1_ftn(message)
785 end function template_fieldname
786 
787 ! ------------------------------------------------------------------------------
788 
789 function field_has_identity(self, fieldname) result(has_identity)
790  class(mpas_geom), intent(in) :: self
791  character(len=*), intent(in) :: fieldname
792  integer :: ii
793  logical :: has_identity
794  has_identity = .false.
795  if (allocated(self%templated_fields)) then
796  do ii = 1, size(self%templated_fields)
797  if (trim(self%templated_fields(ii)%name) == trim(fieldname) .and. &
798  trim(self%templated_fields(ii)%identity) /= 'none') then
799  has_identity = .true.
800  return
801  end if
802  end do
803  end if
804 end function field_has_identity
805 
806 ! ------------------------------------------------------------------------------
807 
808 function identity_fieldname(self, fieldname) result(identity)
809  class(mpas_geom), intent(in) :: self
810  character(len=*), intent(in) :: fieldname
811  integer :: ii
812  character(len=MAXVARLEN) :: identity
813  do ii = 1, size(self%templated_fields)
814  if (trim(self%templated_fields(ii)%name) == trim(fieldname)) then
815  identity = self%templated_fields(ii)%identity
816  return
817  end if
818  end do
819  write(message,*) '--> mpas_geom % identity_fieldname: fieldname does not have identity, ',fieldname
820  call abor1_ftn(message)
821 end function identity_fieldname
822 
823 ! ------------------------------------------------------------------------------
824 
825 !***********************************************************************
826 !
827 ! subroutine pool_has_field
828 !
829 !> \brief Check for presence of field in pool
830 !
831 !-----------------------------------------------------------------------
832 function pool_has_field(pool, field) result(has)
833 
834  implicit none
835 
836  type(mpas_pool_type), pointer, intent(in) :: pool
837  character(len=*), intent(in) :: field
838  type(mpas_pool_iterator_type) :: poolitr
839  logical :: has
840 
841  has = .false.
842 
843  call mpas_pool_begin_iteration(pool)
844 
845  do while ( mpas_pool_get_next_member(pool, poolitr) )
846  if (poolitr % memberType == mpas_pool_field .and. &
847  trim(field) == trim(poolitr % memberName) ) then
848  has = .true.
849  exit
850  endif
851  end do
852 
853 end function pool_has_field
854 
855 
856 !*******************************************************************************
857 !
858 ! function isDecomposed
859 !
860 !> \brief whether horizontal dimension name is decomposed across MPI tasks
861 !
862 !*******************************************************************************
863 
864 logical function isdecomposed(dimName)
865  implicit none
866  character(len=*), intent(in) :: dimname
867  if (trim(dimname) == 'nCells') then
868  isdecomposed = .true.
869  return
870  else if (trim(dimname) == 'nEdges') then
871  isdecomposed = .true.
872  return
873  else if (trim(dimname) == 'nVertices') then
874  isdecomposed = .true.
875  return
876  else
877  isdecomposed = .false.
878  return
879  end if
880 end function isdecomposed
881 
882 
883 !***********************************************************************
884 !
885 ! function getSolveDimNames
886 !
887 !> \brief Returns an array with the dimension names of a pool field member
888 !
889 !-----------------------------------------------------------------------
890 
891 function getsolvedimnames(pool, key) result(solveDimNames)
892 
893  implicit none
894 
895  ! Arguments
896  type(mpas_pool_type), pointer, intent(in) :: pool
897  character(len=*), intent(in) :: key
898 
899  character (len=StrKIND), allocatable :: solvedimnames(:)
900 
901  type(mpas_pool_data_type), pointer :: dptr
902  character (len=StrKIND) :: dimname
903  integer :: id
904 
905  dptr => pool_get_member(pool, key, mpas_pool_field)
906 
907  if (.not. associated(dptr)) then
908  write(message,*)'--> getSolveDimNames: problem finding ',trim(key), &
909  ' MPAS_POOL_FIELD in pool'
910  call abor1_ftn(message)
911  end if
912 
913  allocate(solvedimnames(dptr%contentsDims))
914 
915  do id = 1, dptr%contentsDims
916 
917  ! fields with only one time level
918  !real
919  if (associated(dptr%r1)) then
920  dimname = trim(dptr%r1%dimNames(id))
921  else if (associated(dptr%r2)) then
922  dimname = trim(dptr%r2%dimNames(id))
923  else if (associated(dptr%r3)) then
924  dimname = trim(dptr%r3%dimNames(id))
925  !integer
926  else if (associated(dptr%i1)) then
927  dimname = trim(dptr%i1%dimNames(id))
928  else if (associated(dptr%i2)) then
929  dimname = trim(dptr%i2%dimNames(id))
930  else if (associated(dptr%i3)) then
931  dimname = trim(dptr%i3%dimNames(id))
932 
933  ! fields with multiple time levels (use first time level)
934  !real
935  else if (associated(dptr%r1a)) then
936  dimname = trim(dptr%r1a(1)%dimNames(id))
937  else if (associated(dptr%r2a)) then
938  dimname = trim(dptr%r2a(1)%dimNames(id))
939  else if (associated(dptr%r3a)) then
940  dimname = trim(dptr%r3a(1)%dimNames(id))
941  !integer
942  else if (associated(dptr%i1a)) then
943  dimname = trim(dptr%i1a(1)%dimNames(id))
944  else if (associated(dptr%i2a)) then
945  dimname = trim(dptr%i2a(1)%dimNames(id))
946  else if (associated(dptr%i3a)) then
947  dimname = trim(dptr%i3a(1)%dimNames(id))
948 
949  else
950  write(message,*) '--> getSolveDimNames: unsupported field type for key,&
951  & contentsDims, contentsType = ',trim(key),',',dptr%contentsDims, &
952  ',',dptr%contentsType
953  call abor1_ftn(message)
954  end if
955 
956  if (isdecomposed(dimname)) then
957  solvedimnames(id) = trim(dimname)//'Solve'
958  else
959  solvedimnames(id) = trim(dimname)
960  end if
961  end do
962 end function getsolvedimnames
963 
964 
965 !***********************************************************************
966 !
967 ! function getSolveDimSizes
968 !
969 !> \brief Returns an array with the dimension sizes of a pool field member
970 !
971 !> \details
972 ! In some cases, the pool containing the field will not contain any
973 ! dimensions. Then the optional dimPool argument can be used.
974 !
975 !-----------------------------------------------------------------------
976 
977 function getsolvedimsizes(pool, key, dimPool) result(solveDimSizes)
978 
979  implicit none
980 
981  ! Arguments
982  type (mpas_pool_type), pointer, intent(in) :: pool
983  character(len=*), intent(in) :: key
984  type (mpas_pool_type), pointer, optional, intent(in) :: dimpool
985 
986  integer, allocatable :: solvedimsizes(:)
987 
988  character (len=StrKIND), allocatable :: dimnames(:)
989  integer :: id
990  integer, pointer :: dimsize
991 
992  dimnames = getsolvedimnames(pool, key)
993 
994  allocate(solvedimsizes(size(dimnames)))
995 
996  do id = 1, size(dimnames)
997  if (present(dimpool)) then
998  call mpas_pool_get_dimension(dimpool, dimnames(id), dimsize)
999  else
1000  call mpas_pool_get_dimension(pool, dimnames(id), dimsize)
1001  end if
1002 
1003  if (associated(dimsize)) then
1004  solvedimsizes(id) = dimsize
1005  else
1006  write(message,*)'--> getSolveDimSizes: ',trim(dimnames(id)), &
1007  ' not available'
1008  call abor1_ftn(message)
1009  end if
1010  end do
1011 
1012  deallocate(dimnames)
1013 
1014 end function getsolvedimsizes
1015 
1016 
1017 !***********************************************************************
1018 !
1019 ! function getSolveDimSizes
1020 !
1021 !> \brief Returns the number of vertical levels for a pool field member
1022 !
1023 !> \details
1024 ! In some cases, the pool containing the field will not contain any
1025 ! dimensions. Then the optional dimPool argument can be used.
1026 !
1027 !-----------------------------------------------------------------------
1028 
1029 function getvertlevels(pool, key, dimPool) result(nlevels)
1030 
1031  ! Arguments
1032  type (mpas_pool_type), pointer, intent(in) :: pool
1033  character(len=*), intent(in) :: key
1034  type (mpas_pool_type), pointer, optional, intent(in) :: dimpool
1035 
1036  integer :: nlevels, id
1037  integer, allocatable :: dimsizes(:)
1038  character (len=StrKIND), allocatable :: dimnames(:)
1039 
1040  ! assume default value of 1, unless stated otherwise in field attributes
1041  nlevels = 1
1042 
1043  dimnames = getsolvedimnames(pool, key)
1044  if (present(dimpool)) then
1045  dimsizes = getsolvedimsizes(pool, key, dimpool)
1046  else
1047  dimsizes = getsolvedimsizes(pool, key)
1048  end if
1049  do id = 1, size(dimsizes)
1050  ! check if dimNames(id) is one of the available vertical coordinates
1051  if (ufo_vars_getindex(mpasverticalcoordinates, dimnames(id)) > 0) then
1052  nlevels = dimsizes(id)
1053  end if
1054  end do
1055 
1056  deallocate(dimnames, dimsizes)
1057 
1058 end function getvertlevels
1059 
1060 ! ------------------------------------------------------------------------------
1061 
1062 subroutine var_nlevels(self, var, nlevels)
1063 
1064  implicit none
1065 
1066  class(mpas_geom), intent(in) :: self
1067  character(len=*), intent(in) :: var
1068  integer, intent(out) :: nlevels
1069 
1070  character(len=StrKIND) :: poolVar
1071  type(mpas_pool_type), pointer :: allFields
1072 
1073  allfields => self % domain % blocklist % allFields
1074 
1075  if (pool_has_field(allfields, var)) then
1076  poolvar = var
1077  else if (self % is_templated(var)) then
1078  poolvar = self%template(var)
1079  if (.not.pool_has_field(allfields, poolvar)) then
1080  write(message,*)'--> vars_nlevels: ',trim(poolvar), &
1081  ' not available in MPAS domain'
1082  call abor1_ftn(message)
1083  end if
1084  else
1085  write(message,*)'--> vars_nlevels: ',trim(var), &
1086  ' not available in MPAS domain or self % templated_fields'
1087  call abor1_ftn(message)
1088  end if
1089 
1090  nlevels = getvertlevels(allfields, poolvar, self % domain % blocklist % dimensions)
1091 
1092 end subroutine var_nlevels
1093 
1094 ! ------------------------------------------------------------------------------
1095 
1096 subroutine variables_nlevels(self, vars, nv, nlevels)
1097 
1098  implicit none
1099 
1100  class(mpas_geom), intent(in) :: self
1101  type(oops_variables), intent(in) :: vars
1102  integer(c_size_t), intent(in) :: nv
1103  integer(c_size_t), intent(inout) :: nlevels(nv)
1104 
1105  integer :: ivar, nn
1106 
1107  do ivar = 1, nv
1108  call self%nlevels(vars%variable(ivar), nn)
1109  nlevels(ivar) = int(nn, c_size_t)
1110  end do
1111 
1112 end subroutine variables_nlevels
1113 
1114 ! ------------------------------------------------------------------------------
1115 
1116 end module mpas_geom_mod
real(kind=kind_real), parameter mpas_jedi_half_kr
real(kind=kind_real), parameter mpas_jedi_piio2_kr
real(kind=kind_real), parameter mpas_jedi_rad2deg_kr
subroutine, public geo_clone(self, other)
integer function, dimension(:), allocatable, public getsolvedimsizes(pool, key, dimPool)
Returns an array with the dimension sizes of a pool field member.
subroutine, public geo_delete(self)
subroutine, public geo_info(self, nCellsGlobal, nCells, nCellsSolve, nEdgesGlobal, nEdges, nEdgesSolve, nVertLevels, nVertLevelsP1)
type(registry_t), public mpas_geom_registry
Linked list interface - defines registry_t type.
subroutine, public geo_set_atlas_lonlat(self, afieldset)
character(len=maxvarlen) function template_fieldname(self, fieldname)
character(len=1024) message
subroutine, public geo_is_equal(is_equal, self, other)
type(idcounter), dimension(:), allocatable geom_count
logical function field_is_templated(self, fieldname)
logical function, public pool_has_field(pool, field)
Check for presence of field in pool.
logical function field_has_identity(self, fieldname)
character(len=maxvarlen) function identity_fieldname(self, fieldname)
subroutine var_nlevels(self, var, nlevels)
subroutine, public geo_setup(self, f_conf, f_comm)
Linked list implementation.
character(len=strkind) function, dimension(:), allocatable, public getsolvedimnames(pool, key)
Returns an array with the dimension names of a pool field member.
subroutine variables_nlevels(self, vars, nv, nlevels)
character(len=maxvarlen), dimension(6), parameter mpasverticalcoordinates
logical function isdecomposed(dimName)
whether horizontal dimension name is decomposed across MPI tasks
subroutine, public geo_fill_atlas_fieldset(self, afieldset)
integer function, public getvertlevels(pool, key, dimPool)
Returns the number of vertical levels for a pool field member.
subroutine geo_deallocate_nonda_fields(domain)
Fortran derived type to hold geometry definition.