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
16 use oops_variables_mod,
only: oops_variables
19 use ufo_vars_mod,
only: maxvarlen, ufo_vars_getindex
22 use mpas_derived_types
25 use kinds,
only : kind_real
26 use mpas_dmpar,
only: mpas_dmpar_sum_int
29 use mpas_pool_routines
47 character(len=MAXVARLEN) :: name
48 character(len=MAXVARLEN) :: template
49 character(len=MAXVARLEN) :: identity
54 integer :: ncellsglobal
55 integer :: nedgesglobal
56 integer :: nverticesglobal
60 integer :: ncellssolve
61 integer :: nedgessolve
62 integer :: nverticessolve
63 integer :: nvertlevels
64 integer :: nvertlevelsp1
65 integer :: nsoillevels
66 integer :: vertexdegree
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
86 type (domain_type),
pointer :: domain => null()
87 type (core_type),
pointer :: corelist => null()
89 type(fckit_mpi_comm) :: f_comm
91 type(atlas_functionspace) :: afunctionspace
108 integer :: id, counter
115 character(len=MAXVARLEN),
parameter :: &
117 [character(len=maxvarlen) ::
'nVertLevels',
'nVertLevelsP1', &
118 'nVertLevelsP2',
'nSoilLevels',
'nOznLevels',
'nAerLevels']
120 #define LISTED_TYPE mpas_geom
123 #include <oops/util/linkedList_i.f>
132 #include <oops/util/linkedList_c.f>
140 type(fckit_configuration),
intent(in) :: f_conf
141 type(fckit_mpi_comm),
intent(in) :: f_comm
143 character(len=StrKIND) :: string1
144 type (mpas_pool_type),
pointer :: meshpool, fg
145 type (block_type),
pointer :: block_ptr
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(:)
152 real (kind=kind_real),
pointer :: r1d_ptr(:), r2d_ptr(:,:)
153 integer,
pointer :: i0d_ptr, i1d_ptr(:), i2d_ptr(:,:)
155 type(
idcounter),
allocatable :: prev_count(:)
158 logical :: deallocate_fields
159 logical :: bump_interp
161 call fckit_log%info(
'==> create geom')
167 call f_conf%get_or_die(
"nml_file",str)
169 call f_conf%get_or_die(
"streams_file",str)
173 call mpas_init( self % corelist, self % domain, mpi_comm = self%f_comm%communicator(), &
174 & namelistfileparam = trim(nml_file), streamsfileparam = trim(streams_file))
179 if (f_conf%get(
"interpolation type", str))
then
182 self%use_bump_interpolation = .true.
183 case (
'unstructured')
184 self%use_bump_interpolation = .false.
186 write(
message,*)
'--> geo_setup: interpolation type: ',str,
' not implemented'
190 self%use_bump_interpolation = .true.
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
198 self % deallocate_nonda_fields = .false.
203 if (f_conf%has(
"bump vunit"))
then
204 call f_conf%get_or_die(
"bump vunit",str)
205 self % bump_vunit = str
207 self % bump_vunit =
'modellevel'
212 allocate(prev_count(nprev))
215 if (prev_count(ii)%id == self%domain%domainID)
then
216 call abor1_ftn(
"domainID already used")
231 if (f_conf%get(
'template fields file',str))
then
234 fields_file =
'geovars.yaml'
236 template_conf = fckit_yamlconfiguration(fckit_pathname(fields_file))
237 call template_conf%get_or_die(
'fields',fields_conf)
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)
249 self%templated_fields(ii)%identity =
'none'
252 deallocate(fields_conf)
265 block_ptr => self % domain % blocklist
267 call mpas_pool_get_subpool ( block_ptr % structs,
'mesh', meshpool )
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 )
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 )
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 )
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
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 ) )
321 call mpas_pool_get_array ( meshpool,
'latCell', r1d_ptr )
322 self % latCell = r1d_ptr(1:self % nCells)
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)
365 call mpas_pool_get_array ( meshpool,
'zgrid', r2d_ptr )
366 self % zgrid = r2d_ptr( 1:self % nVertLevelsP1, 1:self % nCells )
368 call fckit_log%debug(
'End of geo_setup')
369 if (
allocated(prev_count))
deallocate(prev_count)
370 if (
allocated(str))
deallocate(str)
381 type(atlas_fieldset),
intent(inout) :: afieldset
383 real(kind_real),
pointer :: real_ptr(:,:)
384 type(atlas_field) :: afield
387 afield = atlas_field(name=
"lonlat", kind=atlas_real(kind_real), shape=(/2,self%nCellsSolve/))
388 call afield%data(real_ptr)
391 call afieldset%add(afield)
402 type(atlas_fieldset),
intent(inout) :: afieldset
405 real(kind=kind_real),
pointer :: real_ptr_1(:), real_ptr_2(:,:)
406 type(atlas_field) :: afield
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)
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
423 ( self%zgrid(jz,1:self%nCellsSolve) + self%zgrid(jz+1,1:self%nCellsSolve) )
430 call afieldset%add(afield)
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
444 type (field0DReal),
pointer :: field0d
445 type (field1DReal),
pointer :: field1d
446 type (field2DReal),
pointer :: field2d
448 integer,
parameter :: num_da_fields = 36
450 character (len=22),
allocatable :: poolname_a(:), poolname_b(:)
451 character (len=22),
allocatable :: da_fieldnames(:)
453 allocate(poolname_a(3))
455 poolname_a(2)=
'tend_physics'
456 poolname_a(3)=
'atm_input'
458 allocate(poolname_b(3))
459 poolname_b(1)=
'diag_physics'
461 poolname_b(3)=
'sfc_input'
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'
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))
511 do i=1,
size(poolname_b)
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)
518 do while ( mpas_pool_get_next_member(pool_b, poolitr_b) )
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
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))
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))
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))
554 deallocate(poolname_a)
555 deallocate(poolname_b)
556 deallocate(da_fieldnames)
572 self%f_comm = other%f_comm
574 call fckit_log%debug(
'====> copy of geom array')
576 self % nCellsGlobal = other % nCellsGlobal
577 self % nCells = other % nCells
578 self % nCellsSolve = other % nCellsSolve
580 self % nEdgesGlobal = other % nEdgesGlobal
581 self % nEdges = other % nEdges
582 self % nEdgesSolve = other % nEdgesSolve
584 self % nVerticesGlobal = other % nVerticesGlobal
585 self % nVertices = other % nVertices
586 self % nVerticesSolve = other % nVerticesSolve
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
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))
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
636 call fckit_log%debug(
'====> copy of geom corelist and domain')
638 self % corelist => other % corelist
639 self % domain => other % domain
640 write(
message,*)
'inside geo_clone: other % domain % domainID = ', other % domain % domainID
644 if (
geom_count(ii)%id == self%domain%domainID)
then
649 call fckit_log%debug(
'====> copy of geom done')
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)
684 if (
geom_count(ii)%id == self%domain%domainID)
then
686 if ((
associated(self % corelist)).and.(
associated(self % domain)))
then
693 '==> destruct MPAS corelist and domain:', self%domain%domainID
695 call mpas_timer_set_context( self % domain )
696 call mpas_finalize(self % corelist, self % domain)
698 nullify(self % corelist)
699 nullify(self % domain)
714 logical*1,
intent(inout) :: is_equal
720 is_equal = (self%nCellsGlobal .eq. other%nCellsGlobal) .and. &
721 (self%nVertLevels .eq. other%nVertLevels)
727 subroutine geo_info(self, nCellsGlobal, nCells, nCellsSolve, &
728 nEdgesGlobal, nEdges, nEdgesSolve, &
729 nVertLevels, nVertLevelsP1)
734 integer,
intent(inout) :: ncellsglobal, ncells, ncellssolve
735 integer,
intent(inout) :: nedgesglobal, nedges, nedgessolve
736 integer,
intent(inout) :: nvertlevels
737 integer,
intent(inout) :: nvertlevelsp1
739 ncellsglobal = self%nCellsGlobal
741 ncellssolve = self%nCellsSolve
743 nedgesglobal = self%nEdgesGlobal
745 nedgessolve = self%nEdgesSolve
747 nvertlevels = self%nVertLevels
748 nvertlevelsp1 = self%nVertLevelsP1
756 character(len=*),
intent(in) :: fieldname
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.
774 character(len=*),
intent(in) :: fieldname
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
783 write(
message,*)
'--> mpas_geom % template_fieldname: fieldname is not templated, ',fieldname
791 character(len=*),
intent(in) :: fieldname
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.
810 character(len=*),
intent(in) :: fieldname
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
819 write(
message,*)
'--> mpas_geom % identity_fieldname: fieldname does not have identity, ',fieldname
836 type(mpas_pool_type),
pointer,
intent(in) :: pool
837 character(len=*),
intent(in) :: field
838 type(mpas_pool_iterator_type) :: poolitr
843 call mpas_pool_begin_iteration(pool)
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
866 character(len=*),
intent(in) :: dimname
867 if (trim(dimname) ==
'nCells')
then
870 else if (trim(dimname) ==
'nEdges')
then
873 else if (trim(dimname) ==
'nVertices')
then
896 type(mpas_pool_type),
pointer,
intent(in) :: pool
897 character(len=*),
intent(in) :: key
899 character (len=StrKIND),
allocatable :: solvedimnames(:)
901 type(mpas_pool_data_type),
pointer :: dptr
902 character (len=StrKIND) :: dimname
905 dptr => pool_get_member(pool, key, mpas_pool_field)
907 if (.not.
associated(dptr))
then
908 write(
message,*)
'--> getSolveDimNames: problem finding ',trim(key), &
909 ' MPAS_POOL_FIELD in pool'
913 allocate(solvedimnames(dptr%contentsDims))
915 do id = 1, dptr%contentsDims
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))
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))
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))
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))
950 write(
message,*)
'--> getSolveDimNames: unsupported field type for key,&
951 & contentsDims, contentsType = ',trim(key),
',',dptr%contentsDims, &
952 ',',dptr%contentsType
957 solvedimnames(id) = trim(dimname)//
'Solve'
959 solvedimnames(id) = trim(dimname)
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
986 integer,
allocatable :: solvedimsizes(:)
988 character (len=StrKIND),
allocatable :: dimnames(:)
990 integer,
pointer :: dimsize
994 allocate(solvedimsizes(
size(dimnames)))
996 do id = 1,
size(dimnames)
997 if (
present(dimpool))
then
998 call mpas_pool_get_dimension(dimpool, dimnames(id), dimsize)
1000 call mpas_pool_get_dimension(pool, dimnames(id), dimsize)
1003 if (
associated(dimsize))
then
1004 solvedimsizes(id) = dimsize
1006 write(
message,*)
'--> getSolveDimSizes: ',trim(dimnames(id)), &
1012 deallocate(dimnames)
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
1036 integer :: nlevels, id
1037 integer,
allocatable :: dimsizes(:)
1038 character (len=StrKIND),
allocatable :: dimnames(:)
1044 if (
present(dimpool))
then
1049 do id = 1,
size(dimsizes)
1052 nlevels = dimsizes(id)
1056 deallocate(dimnames, dimsizes)
1067 character(len=*),
intent(in) :: var
1068 integer,
intent(out) :: nlevels
1070 character(len=StrKIND) :: poolVar
1071 type(mpas_pool_type),
pointer :: allFields
1073 allfields => self % domain % blocklist % allFields
1077 else if (self % is_templated(var))
then
1078 poolvar = self%template(var)
1080 write(
message,*)
'--> vars_nlevels: ',trim(poolvar), &
1081 ' not available in MPAS domain'
1085 write(
message,*)
'--> vars_nlevels: ',trim(var), &
1086 ' not available in MPAS domain or self % templated_fields'
1090 nlevels =
getvertlevels(allfields, poolvar, self % domain % blocklist % dimensions)
1101 type(oops_variables),
intent(in) :: vars
1102 integer(c_size_t),
intent(in) :: nv
1103 integer(c_size_t),
intent(inout) :: nlevels(nv)
1108 call self%nlevels(vars%variable(ivar), nn)
1109 nlevels(ivar) = int(nn, c_size_t)
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.