MPAS-JEDI
mpas_increment_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_geometry, atlas_field, atlas_fieldset, atlas_real
9 use fckit_configuration_module, only: fckit_configuration
10 use fckit_log_module, only: fckit_log
11 
12 !oops
13 use datetime_mod
14 use kinds, only: kind_real
15 use oops_variables_mod
16 
17 !ufo
18 use ufo_geovals_mod
19 use ufo_vars_mod
20 
21 !MPAS-Model
22 use mpas_constants
23 use mpas_dmpar
24 use mpas_derived_types
25 use mpas_kind_types, only: strkind
26 use mpas_pool_routines
27 
28 !mpas-jedi
30 use mpas_geom_mod
33 use mpas4da_mod
34 
35 implicit none
36 
37 private
38 
39 public :: diff_incr, dirac, &
41 
42 integer, parameter :: max_string=8000
43 character(max_string) :: message
44 
45 ! ------------------------------------------------------------------------------
46 
47 contains
48 
49 ! ------------------------------------------------------------------------------
50 
51 subroutine diff_incr(lhs,x1,x2)
52 
53  implicit none
54  class(mpas_fields), intent(inout) :: lhs
55  class(mpas_fields), intent(in) :: x1
56  class(mpas_fields), intent(in) :: x2
57  character(len=StrKIND) :: kind_op
58 
59  call lhs%zeros()
60  if (x1%geom%nCells==x2%geom%nCells .and. x1%geom%nVertLevels==x2%geom%nVertLevels) then
61  if (lhs%geom%nCells==x1%geom%nCells .and. lhs%geom%nVertLevels==x1%geom%nVertLevels) then
62  kind_op = 'sub'
63  call da_operator(trim(kind_op), lhs % subFields, x1 % subFields, x2 % subFields, fld_select = lhs % fldnames_ci)
64  else
65  call abor1_ftn("mpas_increment:diff_incr: dimension mismatch between the two variables.")
66  endif
67  else
68  call abor1_ftn("mpas_increment:diff_incr: states not at same resolution")
69  endif
70 
71  return
72 
73 end subroutine diff_incr
74 
75 ! ------------------------------------------------------------------------------
76 
77 subroutine dirac(self, f_conf)
78 
79  implicit none
80  class(mpas_fields), intent(inout) :: self
81  type(fckit_configuration), intent(in) :: f_conf !< Configuration
82  character(len=:), allocatable :: str
83  integer :: ndir, idir, ildir, ndirlocal
84  character(len=3) :: idirchar
85  character(len=StrKIND) :: dirvar
86  type (mpas_pool_iterator_type) :: poolitr
87  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a
88  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a
89  integer :: nearestcell
90  integer, allocatable, dimension(:) :: dirowned, dirownedglobal
91  real (kind=kind_real), allocatable, dimension(:) :: dirlats
92  real (kind=kind_real), allocatable, dimension(:) :: dirlons
93  integer, allocatable, dimension(:) :: dircells
94 
95  ! Get number and positions of Diracs
96  call f_conf%get_or_die("ndir",ndir)
97 
98  allocate( dirowned(ndir) )
99  allocate( dircells(ndir) )
100 
101  call f_conf%get_or_die("dirLats",dirlats)
102  call f_conf%get_or_die("dirLons",dirlons)
103 
104  call f_conf%get_or_die("ildir",ildir)
105  call f_conf%get_or_die("dirvar",str)
106  dirvar = str
107 
108  !Test if dir is owned and find the nearest local cell
109  ! (repurposed from MPAS-Release/src/core_atmosphere/diagnostics/soundings.F)
110 
111  ndirlocal = 0
112  do idir=1,ndir
113  nearestcell = self % geom % nCellsSolve
114  nearestcell = nearest_cell( (dirlats(idir) * mpas_jedi_deg2rad_kr), &
115  (dirlons(idir) * mpas_jedi_deg2rad_kr), &
116  nearestcell, self % geom % nCells, self % geom % maxEdges, &
117  self % geom % nEdgesOnCell, self % geom % cellsOnCell, &
118  self % geom % latCell, self % geom % lonCell )
119 
120  if (nearestcell <= self % geom % nCellsSolve) then
121  dirowned(idir) = 1
122  dircells(idir) = nearestcell
123  ndirlocal = ndirlocal + 1
124  else
125  dirowned(idir) = 0
126  dircells(idir) = self % geom % nCells + 1
127  end if
128  end do
129 
130  write(message,*) 'This processor owns ',ndirlocal,' dirac forcing locations'
131  call fckit_log%debug(message)
132 
133  ! Check
134  if (ndir<1) call abor1_ftn("mpas_increment:dirac non-positive ndir")
135 
136  allocate( dirownedglobal(ndir) )
137  call mpas_dmpar_max_int_array( self % geom % domain % dminfo, ndir, dirowned, dirownedglobal)
138  if ( any(dirownedglobal.lt.1) ) then
139  call abor1_ftn("mpas_increment:dirac invalid Lat/Lon")
140  end if
141 
142  call mpas_dmpar_sum_int_array( self % geom % domain % dminfo, ndir, dirowned, dirownedglobal)
143  if ( any(dirownedglobal.gt.1) ) then
144  call abor1_ftn("mpas_increment:duplicated dirac on >1 processors")
145  end if
146  deallocate( dirownedglobal )
147 
148  if ((ildir < 1) .or. (ildir > self % geom % nVertLevels)) then
149  call abor1_ftn("mpas_increment:dirac invalid ildir")
150  endif
151 
152  ! Setup Diracs
153  call self%zeros()
154 
155  call mpas_pool_begin_iteration(self % subFields)
156 
157  do while ( mpas_pool_get_next_member(self % subFields, poolitr) )
158  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
159  ! so we select only those members of the pool that are fields
160  if (poolitr % memberType == mpas_pool_field) then
161  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
162  write(message,*) 'poolItr % nDims , poolItr % memberName =', poolitr % nDims , poolitr % memberName
163  call fckit_log%debug(message)
164 
165  if (poolitr % dataType == mpas_pool_real) then
166  ! Depending on the dimensionality of the field, we need to set pointers of
167  ! the correct type
168  if( trim(dirvar) .eq. trim(poolitr % memberName) ) then
169  ndirlocal = 0
170  do idir=1, ndir
171  if ( dirowned(idir).eq.1 ) then
172  if (poolitr % nDims == 1) then
173  call mpas_pool_get_array(self % subFields, trim(poolitr % memberName), r1d_ptr_a)
174  r1d_ptr_a( dircells(idir) ) = mpas_jedi_one_kr
175  else if (poolitr % nDims == 2) then
176  call mpas_pool_get_array(self % subFields, trim(poolitr % memberName), r2d_ptr_a)
177  r2d_ptr_a( ildir, dircells(idir) ) = mpas_jedi_one_kr
178  else if (poolitr % nDims == 3) then
179  call fckit_log%info ('Not implemented yet')
180  end if
181  ndirlocal = ndirlocal + 1
182  end if
183  end do
184  end if
185  end if
186  end if
187  end do
188 
189  deallocate( dirowned )
190  deallocate( dirlats )
191  deallocate( dirlons )
192  deallocate( dircells )
193 
194 end subroutine dirac
195 
196 ! ------------------------------------------------------------------------------
197 
198 !!TODO: Alternatively could make the function nearest_cell public in MPAS
199 !! and then define an interface to it within this module (cleaner)
200 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
201 ! Finds the MPAS grid cell nearest to (target_lat, target_lon)
202 !!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
203 integer function nearest_cell(target_lat, target_lon, start_cell, nCells, maxEdges, &
204  nEdgesOnCell, cellsOnCell, latCell, lonCell)
205 
206  implicit none
207 
208  real (kind=kind_real), intent(in) :: target_lat, target_lon
209  integer, intent(in) :: start_cell
210  integer, intent(in) :: ncells, maxedges
211  integer, dimension(nCells), intent(in) :: nedgesoncell
212  integer, dimension(maxEdges,nCells), intent(in) :: cellsoncell
213  real (kind=kind_real), dimension(nCells), intent(in) :: latcell, loncell
214 
215  integer :: i
216  integer :: icell
217  integer :: current_cell
218  real (kind=kind_real) :: current_distance, d
219  real (kind=kind_real) :: nearest_distance
220  type (atlas_geometry) :: ageometry
221 
222  nearest_cell = start_cell
223  current_cell = -1
224 
225  ageometry = atlas_geometry("UnitSphere")
226 
227  do while (nearest_cell /= current_cell)
228  current_cell = nearest_cell
229  current_distance = ageometry%distance(loncell(current_cell), latcell(current_cell), target_lon, target_lat)
230  nearest_cell = current_cell
231  nearest_distance = current_distance
232  do i = 1, nedgesoncell(current_cell)
233  icell = cellsoncell(i,current_cell)
234  if (icell <= ncells) then
235  d = ageometry%distance(loncell(icell), latcell(icell), target_lon, target_lat)
236  if (d < nearest_distance) then
237  nearest_cell = icell
238  nearest_distance = d
239  end if
240  end if
241  end do
242  end do
243 end function nearest_cell
244 
245 ! ------------------------------------------------------------------------------
246 
247 subroutine set_atlas(self, geom, vars, afieldset)
248 
249  implicit none
250 
251  type(mpas_fields), intent(in) :: self
252  type(mpas_geom), intent(in) :: geom
253  type(oops_variables), intent(in) :: vars
254  type(atlas_fieldset), intent(inout) :: afieldset
255 
256  integer :: jvar, nlevels
257  logical :: var_found
258  type(atlas_field) :: afield
259  type(mpas_pool_iterator_type) :: poolitr
260 
261  do jvar = 1,vars%nvars()
262  var_found = .false.
263  call mpas_pool_begin_iteration(self%subFields)
264  do while (mpas_pool_get_next_member(self%subFields, poolitr))
265  if (trim(vars%variable(jvar))==trim(poolitr%memberName)) then
266  if (.not.afieldset%has_field(vars%variable(jvar))) then
267  ! Create field
268  if (poolitr%nDims==1) then
269  nlevels = 0
270  else if (poolitr%nDims==2) then
271  nlevels = geom%nVertLevels
272  else if (poolitr%nDims==3) then
273  call abor1_ftn('not implemented yet')
274  end if
275  afield = geom%afunctionspace%create_field(name=vars%variable(jvar),kind=atlas_real(kind_real),levels=nlevels)
276 
277  ! Add field
278  call afieldset%add(afield)
279 
280  ! Release pointer
281  call afield%final()
282  end if
283 
284  ! Set flag
285  var_found = .true.
286  exit
287  end if
288  end do
289  if (.not.var_found) call abor1_ftn('variable '//trim(vars%variable(jvar))//' not found in increment')
290  end do
291 
292 end subroutine set_atlas
293 
294 ! ------------------------------------------------------------------------------
295 
296 subroutine to_atlas(self, geom, vars, afieldset)
297 
298  implicit none
299 
300  type(mpas_fields), intent(in) :: self
301  type(mpas_geom), intent(in) :: geom
302  type(oops_variables), intent(in) :: vars
303  type(atlas_fieldset), intent(inout) :: afieldset
304 
305  integer :: jvar, nlevels
306  real(kind=kind_real), pointer :: real_ptr_1(:), real_ptr_2(:,:)
307  real(kind=kind_real), pointer :: r1d_ptr_a(:), r2d_ptr_a(:,:)
308  logical :: var_found
309  type(atlas_field) :: afield
310  type(mpas_pool_iterator_type) :: poolitr
311 
312  do jvar = 1,vars%nvars()
313  var_found = .false.
314  call mpas_pool_begin_iteration(self%subFields)
315  do while (mpas_pool_get_next_member(self%subFields, poolitr))
316  if (trim(vars%variable(jvar))==trim(poolitr%memberName)) then
317  if (afieldset%has_field(vars%variable(jvar))) then
318  ! Get field
319  afield = afieldset%field(vars%variable(jvar))
320  else
321  ! Create field
322  if (poolitr%nDims==1) then
323  nlevels = 0
324  else if (poolitr%nDims==2) then
325  nlevels = geom%nVertLevels
326  else if (poolitr%nDims==3) then
327  call abor1_ftn('not implemented yet')
328  end if
329  afield = geom%afunctionspace%create_field(name=vars%variable(jvar),kind=atlas_real(kind_real),levels=nlevels)
330 
331  ! Add field
332  call afieldset%add(afield)
333  end if
334 
335  ! Copy data
336  if (poolitr%nDims==1) then
337  call afield%data(real_ptr_1)
338  call mpas_pool_get_array(self%subFields, trim(poolitr%memberName), r1d_ptr_a)
339  real_ptr_1 = r1d_ptr_a(1:geom%nCellsSolve)
340  else if (poolitr%nDims==2) then
341  call afield%data(real_ptr_2)
342  call mpas_pool_get_array(self%subFields, trim(poolitr%memberName), r2d_ptr_a)
343  real_ptr_2 = r2d_ptr_a(1:geom%nVertLevels,1:geom%nCellsSolve)
344  end if
345 
346  ! Release pointer
347  call afield%final()
348 
349  ! Set flag
350  var_found = .true.
351  exit
352  end if
353  end do
354  if (.not.var_found) call abor1_ftn('variable '//trim(vars%variable(jvar))//' not found in increment')
355  end do
356 
357 end subroutine to_atlas
358 
359 ! ------------------------------------------------------------------------------
360 
361 subroutine from_atlas(self, geom, vars, afieldset)
362 
363  implicit none
364 
365  type(mpas_fields), intent(inout) :: self
366  type(mpas_geom), intent(in) :: geom
367  type(oops_variables), intent(in) :: vars
368  type(atlas_fieldset), intent(in) :: afieldset
369 
370  integer :: jvar
371  real(kind=kind_real), pointer :: real_ptr_1(:), real_ptr_2(:,:)
372  real(kind=kind_real), pointer :: r1d_ptr_a(:), r2d_ptr_a(:,:)
373  logical :: var_found
374  type(atlas_field) :: afield
375  type(mpas_pool_iterator_type) :: poolitr
376 
377  do jvar = 1,vars%nvars()
378  var_found = .false.
379  call mpas_pool_begin_iteration(self%subFields)
380  do while (mpas_pool_get_next_member(self%subFields, poolitr))
381  if (trim(vars%variable(jvar))==trim(poolitr%memberName)) then
382  ! Get field
383  afield = afieldset%field(vars%variable(jvar))
384 
385  ! Copy data
386  if (poolitr%nDims==1) then
387  call afield%data(real_ptr_1)
388  call mpas_pool_get_array(self%subFields, trim(poolitr%memberName), r1d_ptr_a)
389  r1d_ptr_a(1:geom%nCellsSolve) = real_ptr_1
390  else if (poolitr%nDims==2) then
391  call afield%data(real_ptr_2)
392  call mpas_pool_get_array(self%subFields, trim(poolitr%memberName), r2d_ptr_a)
393  r2d_ptr_a(1:geom%nVertLevels,1:geom%nCellsSolve) = real_ptr_2
394  else if (poolitr%nDims==3) then
395  call abor1_ftn('not implemented yet')
396  end if
397 
398  ! Release pointer
399  call afield%final()
400 
401  ! Set flag
402  var_found = .true.
403  exit
404  end if
405  end do
406  if (.not.var_found) call abor1_ftn('variable '//trim(vars%variable(jvar))//' not found in increment')
407  end do
408 
409  ! TODO: Since only local locations are updated/transferred from ug,
410  ! need MPAS HALO comms before using these fields in MPAS
411 
412 end subroutine from_atlas
413 
414 ! ------------------------------------------------------------------------------
415 
416 end module mpas_increment_mod
integer, parameter max_string
character(max_string) message
subroutine, public da_operator(kind_op, pool_a, pool_b, pool_c, fld_select)
Performs A = A 'kind_op' B for pools A and B.
real(kind=kind_real), parameter mpas_jedi_deg2rad_kr
real(kind=kind_real), parameter mpas_jedi_one_kr
subroutine, public set_atlas(self, geom, vars, afieldset)
subroutine, public to_atlas(self, geom, vars, afieldset)
integer function nearest_cell(target_lat, target_lon, start_cell, nCells, maxEdges, nEdgesOnCell, cellsOnCell, latCell, lonCell)
subroutine, public diff_incr(lhs, x1, x2)
subroutine, public dirac(self, f_conf)
subroutine, public from_atlas(self, geom, vars, afieldset)
Fortran derived type to hold MPAS field.
Fortran derived type to hold geometry definition.