Loading [MathJax]/jax/input/TeX/config.js
MPAS-JEDI
All Classes Namespaces Files Functions Variables Typedefs Macros Pages
mpas4da_mod.F90
Go to the documentation of this file.
1 ! Copyright (c) 2018, National Atmospheric for Atmospheric Research (NCAR).
2 !
3 ! Unless noted otherwise source code is licensed under the BSD license.
4 ! Additional copyright and license information can be found in the LICENSE file
5 ! distributed with this code, or at http://mpas-dev.github.com/license.html
6 
7 module mpas4da_mod
8 
9  !***********************************************************************
10  !
11  ! Module mpas4da_mod to encapsulate operations needed for
12  ! Data assimilation purpose.
13  ! It can be used from /somewhere/MPAS/src/operators
14  ! or from /somewhere/mpas-bundle/mpas/model (OOPS)
15  !> \author Gael Descombes/Mickael Duda NCAR/MMM
16  !> \date January 2018
17  !
18  !-----------------------------------------------------------------------
19 
20 !fckit
21 use fckit_log_module, only: fckit_log
22 
23 !oops
24 use kinds, only: kind_real
25 
26 !saber?
27 use random_mod, only: normal_distribution
28 
29 !ufo
30 use ufo_vars_mod
31 
32 !MPAS-Model
33 use mpas_abort, only: mpas_dmpar_global_abort
34 use mpas_constants
35 use mpas_derived_types
36 use mpas_dmpar
37 use mpas_field_routines
38 use mpas_kind_types, only: shortstrkind
39 use mpas_pool_routines
40 
41 !mpas-jedi
44 
45 private
46 
47 public :: &
52  !mpas_pool_template_field, &
53  da_random, &
54  da_operator, &
55  da_self_mult, &
56  da_constant, &
57  da_posdef, &
58  da_setval, &
59  da_axpy, &
60  da_gpnorm, &
61  da_fldrms, &
66 
67 character(len=1024) :: message
68 
69 contains
70 
71  !***********************************************************************
72  !
73  ! function field_is_scalar
74  !
75  !> \brief Test for a form of water vapor or hydrometeor of interest
76  !> \details
77  !> At various places in this module we wish to test for the case where
78  !> a string is one of several 'q?' values. Rather than repeat that
79  !> logic many times, it is encapsulated here.
80  !
81  !-----------------------------------------------------------------------
82  pure function field_is_scalar(fieldName)
83 
84  implicit none
85 
86  character (len=*), intent(in) :: fieldname
87  logical :: field_is_scalar
88  field_is_scalar = any(trim(fieldname) == &
89  (/'qv', 'qc', 'qi', 'qr', 'qs', 'qg', 'qh', 'nr'/))
90 
91  end function
92 
93  !***********************************************************************
94  !
95  ! function match_scalar
96  !
97  !> \brief Test for a form of water vapor or hydrometeor of interest
98  !> \author Steven Vahl
99  !> \date 11 July 2019
100  !> \details
101  !> Test for the case where one string is 'scalars' and another string
102  !> matches with available scalar variables.
103  !
104  !-----------------------------------------------------------------------
105  pure function match_scalar(scalarName, fieldName)
106 
107  implicit none
108 
109  character (len=*), intent(in) :: scalarname
110  character (len=*), intent(in) :: fieldname
111  logical :: match_scalar
112 
113  match_scalar = ( &
114  scalarname == 'scalars' .and. &
115  field_is_scalar(fieldname) )
116 
117  end function
118 
119  !***********************************************************************
120  !
121  ! subroutine mpas_pool_demo
122  !
123  !> \brief Demonstrate basic usage of MPAS pools
124  !> \author Michael Duda
125  !> \date 20 December 2017
126  !> \details
127  !> This routine provides a simple demonstration of how to construct a new
128  !> pool at runtime, add members (fields) to the pool, and to perform generic
129  !> operations on that pool.
130  !
131  !-----------------------------------------------------------------------
132  subroutine mpas_pool_demo(block)
133 
134  implicit none
135 
136  type (block_type), pointer :: block
137 
138  type (mpas_pool_type), pointer :: structs
139  type (mpas_pool_type), pointer :: allFields
140  type (mpas_pool_type), pointer :: da_state
141  type (mpas_pool_type), pointer :: da_state_incr
142 
143  type (field2DReal), pointer :: field
144 
145  write(0,*) '****** Begin pool demo routine ******'
146 
147  structs => block % structs
148  allfields => block % allFields
149 
150  !
151  ! Create a new pool
152  !
153  call mpas_pool_create_pool(da_state)
154 
155  !
156  ! Get pointers to several fields from the allFields pool, and add
157  ! those fields to the da_state pool as well
158  !
159  call mpas_pool_get_field(allfields, 'theta', field)
160  call mpas_pool_add_field(da_state, 'theta', field)
161  write(0,*) 'Now, max value of theta is ', maxval(field % array),minval(field % array)
162  field % array(:,:) = mpas_jedi_one_kr
163  write(0,*)'Dimensions Field: ',field % dimSizes(:)
164 
165  call mpas_pool_get_field(allfields, 'rho', field)
166  call mpas_pool_add_field(da_state, 'rho', field)
167  write(0,*) 'Now, max value of rho is ', maxval(field % array),minval(field % array)
168  field % array(:,:) = mpas_jedi_one_kr
169 
170  !
171  ! Create another pool
172  !
173  call mpas_pool_create_pool(da_state_incr)
174 
175  !
176  ! Duplicate the members of da_state into da_state_incr, and do a deep
177  ! copy of the fields from da_state to da_state_incr
178  !
179  call mpas_pool_clone_pool(da_state, da_state_incr)
180 
181  !
182  ! Call example algebra routine to compute A = A + B for all fields in
183  ! the da_state and da_state_inc pools
184  !
185  call da_operator_addition(da_state, da_state_incr)
186 
187  call mpas_pool_get_field(da_state_incr, 'rho', field)
188  write(0,*) 'Now, max value of rho_incr is ', maxval(field % array)
189 
190  call mpas_pool_get_field(da_state, 'rho', field)
191  write(0,*) 'Now, max value of rho is ', maxval(field % array)
192 
193  !
194  ! Before destroying a pool, we should remove any fields that are
195  ! still referenced by other active pools to avoid deallocating them
196  !
197  call mpas_pool_empty_pool(da_state)
198 
199  !
200  ! Destroy the now-empty da_state pool
201  !
202  call mpas_pool_destroy_pool(da_state)
203 
204  !
205  ! Destroy the da_state_incr pool, deallocating all of its
206  ! fields in the process (because this pool was not emptied)
207  !
208  call mpas_pool_destroy_pool(da_state_incr)
209 
210  write(0,*) '****** End pool demo routine ******'
211 
212  end subroutine mpas_pool_demo
213 
214  !***********************************************************************
215  !
216  ! subroutine da_operator_addition
217  !
218  !> \brief Performs A = A + B for pools A and B
219  !> \author Michael Duda
220  !> \date 20 December 2017
221  !> \details
222  !> Given two pools, A and B, where the fields in B are a subset of
223  !> the fields in A, this routine adds the fields in B to fields in A
224  !> with the same name. When A and B contain identical fields, this
225  !> is equivalent to A = A + B.
226  !
227  !-----------------------------------------------------------------------
228  subroutine da_operator_addition(pool_a, pool_b)
229 
230  implicit none
231 
232  type (mpas_pool_type), pointer :: pool_a, pool_b
233 
234  type (mpas_pool_iterator_type) :: poolitr
235  real (kind=kind_real), pointer :: r0d_ptr_a, r0d_ptr_b
236  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a, r1d_ptr_b
237  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a, r2d_ptr_b
238  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a, r3d_ptr_b
239 
240  !
241  ! Iterate over all fields in pool_b, adding them to fields of the same
242  ! name in pool_a
243  !
244  call mpas_pool_begin_iteration(pool_b)
245 
246  do while ( mpas_pool_get_next_member(pool_b, poolitr) )
247 
248  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
249  ! so we select only those members of the pool that are fields
250  if (poolitr % memberType == mpas_pool_field) then
251 
252  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
253  if (poolitr % dataType == mpas_pool_real) then
254 
255  ! Depending on the dimensionality of the field, we need to set pointers of
256  ! the correct type
257  if (poolitr % nDims == 0) then
258  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
259  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r0d_ptr_b)
260  r0d_ptr_a = r0d_ptr_a + r0d_ptr_b
261  else if (poolitr % nDims == 1) then
262  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
263  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r1d_ptr_b)
264  r1d_ptr_a = r1d_ptr_a + r1d_ptr_b
265  else if (poolitr % nDims == 2) then
266  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
267  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r2d_ptr_b)
268  r2d_ptr_a = r2d_ptr_a + r2d_ptr_b
269  write(message,*) 'Operator add MIN/MAX: ',minval(r2d_ptr_a),maxval(r2d_ptr_a)
270  call fckit_log%debug(message)
271  else if (poolitr % nDims == 3) then
272  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
273  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r3d_ptr_b)
274  r3d_ptr_a = r3d_ptr_a + r3d_ptr_b
275  end if
276 
277  end if
278  end if
279  end do
280 
281  end subroutine da_operator_addition
282 
283 
284  !***********************************************************************
285  !
286  ! subroutine da_copy_all2sub_fields
287  !
288  !> \brief Performs a copy of allfield to a sub pool A
289  !> \author Gael Desccombes
290  !> \date 5 February 2018
291  !> \details
292  !> Given two pools, allfields and A, where the fields in A are a subset of
293  !> the fields in allfields, this routine copy the fields allfields to fields in A
294  !> with the same name.
295  !
296  !-----------------------------------------------------------------------
297  subroutine da_copy_all2sub_fields(domain, pool_a)
298 
299  implicit none
300 
301  type (mpas_pool_type), pointer, intent(inout) :: pool_a
302  type (mpas_pool_type), pointer :: pool_b, state
303  type (domain_type), pointer, intent(in) :: domain
304 
305  type (mpas_pool_iterator_type) :: poolitr_a, poolitr_b
306  real (kind=kind_real), pointer :: r0d_ptr_a, r0d_ptr_b
307  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a, r1d_ptr_b
308  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a, r2d_ptr_b
309  integer, pointer :: index_scalar
310 
311  type (field2dreal), pointer :: field2d
312  type (field3dreal), pointer :: field3d
313 
314 
315  pool_b => domain % blocklist % allFields
316  call mpas_pool_get_subpool(domain % blocklist % structs,'state',state)
317  !
318  ! Iterate over all fields in pool_b, adding them to fields of the same
319  ! name in pool_a
320  !
321  call mpas_pool_begin_iteration(pool_b)
322 
323  do while ( mpas_pool_get_next_member(pool_b, poolitr_b) )
324 
325  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
326  ! so we select only those members of the pool that are fields
327  if (poolitr_b % memberType == mpas_pool_field) then
328 
329  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
330  if (poolitr_b % dataType == mpas_pool_real) then
331 
332  call mpas_pool_begin_iteration(pool_a)
333  do while ( mpas_pool_get_next_member(pool_a, poolitr_a) )
334 
335  if (( trim(poolitr_b % memberName)).eq.(trim(poolitr_a % memberName)) ) then
336 
337  ! Depending on the dimensionality of the field, we need to set pointers of
338  ! the correct type
339  if (poolitr_b % nDims == 0) then
340  call mpas_pool_get_array(pool_a, trim(poolitr_a % memberName), r0d_ptr_a)
341  call mpas_pool_get_array(pool_b, trim(poolitr_b % memberName), r0d_ptr_b)
342  r0d_ptr_a = r0d_ptr_b
343  else if (poolitr_b % nDims == 1) then
344  call mpas_pool_get_array(pool_a, trim(poolitr_a % memberName), r1d_ptr_a)
345  call mpas_pool_get_array(pool_b, trim(poolitr_b % memberName), r1d_ptr_b)
346  r1d_ptr_a = r1d_ptr_b
347  else if (poolitr_b % nDims == 2) then
348  write(message,*) 'poolItr_b % memberName=',trim(poolitr_b % memberName)
349  call fckit_log%debug(message)
350  call mpas_pool_get_array(pool_a, trim(poolitr_a % memberName), r2d_ptr_a)
351  call mpas_pool_get_array(pool_b, trim(poolitr_b % memberName), r2d_ptr_b)
352  r2d_ptr_a = r2d_ptr_b
353  write(message,*) 'Copy all2sub field MIN/MAX: ',trim(poolitr_b % memberName), &
354  minval(r2d_ptr_a),maxval(r2d_ptr_a)
355  call fckit_log%debug(message)
356  end if
357 
358  else if ( match_scalar(trim(poolitr_b % memberName), trim(poolitr_a % memberName)) ) then
359  write(message,*) 'Copy all2sub field: Looking at SCALARS now',trim(poolitr_a % memberName)
360  call fckit_log%debug(message)
361  call mpas_pool_get_dimension(state, 'index_'//trim(poolitr_a % memberName), index_scalar)
362  if (index_scalar .gt. 0) then
363  call mpas_pool_get_field(pool_a, trim(poolitr_a % memberName), field2d)
364  call mpas_pool_get_field(pool_b, trim(poolitr_b % memberName), field3d)
365  field2d % array(:,:) = field3d % array(index_scalar,:,:)
366  write(message,*) 'Copy all2sub field MIN/MAX: ',trim(poolitr_a % memberName), &
367  minval(field2d % array), maxval(field2d % array)
368  call fckit_log%debug(message)
369  else
370  write(message,*) 'WARNING in Copy all2sub field; ',trim(poolitr_a % memberName), &
371  'not available from MPAS'
372  call fckit_log%debug(message)
373  end if
374  end if
375  end do
376  end if
377  end if
378  end do
379 
380  end subroutine da_copy_all2sub_fields
381 
382 
383  !***********************************************************************
384  !
385  ! subroutine da_copy_sub2all_fields
386  !
387  !> \brief Performs a copy of a sub pool A to allfields
388  !> \author Gael Desccombes
389  !> \date 5 February 2018
390  !> \details
391  !> Given two pools, allfields and A, where the fields in A are a subset of
392  !> the fields in allfields, this routine copy the subfields to allfields
393  !> with the same name.
394  !
395  !-----------------------------------------------------------------------
396  subroutine da_copy_sub2all_fields(domain, pool_a)
397 
398  implicit none
399 
400  type (mpas_pool_type), pointer, intent(in) :: pool_a
401  type (mpas_pool_type), pointer :: pool_b, state
402  type (domain_type), pointer, intent(inout) :: domain
403 
404  type (mpas_pool_iterator_type) :: poolitr_a, poolitr_b
405  real (kind=kind_real), pointer :: r0d_ptr_a, r0d_ptr_b
406  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a, r1d_ptr_b
407  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a, r2d_ptr_b
408  integer, pointer :: index_scalar
409 
410  type (field2dreal), pointer :: field2d
411  type (field3dreal), pointer :: field3d
412 
413 
414  pool_b => domain % blocklist % allFields
415  call mpas_pool_get_subpool(domain % blocklist % structs,'state',state)
416  !
417  ! Iterate over all fields in pool_b, adding them to fields of the same
418  ! name in pool_a
419  !
420  call mpas_pool_begin_iteration(pool_b)
421 
422  do while ( mpas_pool_get_next_member(pool_b, poolitr_b) )
423 
424  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
425  ! so we select only those members of the pool that are fields
426  if (poolitr_b % memberType == mpas_pool_field) then
427 
428  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
429  if (poolitr_b % dataType == mpas_pool_real) then
430 
431  call mpas_pool_begin_iteration(pool_a)
432  do while ( mpas_pool_get_next_member(pool_a, poolitr_a) )
433 
434  if (( trim(poolitr_b % memberName)).eq.(trim(poolitr_a % memberName)) ) then
435 
436  ! Depending on the dimensionality of the field, we need to set pointers of
437  ! the correct type
438  if (poolitr_b % nDims == 0) then
439  call mpas_pool_get_array(pool_a, trim(poolitr_a % memberName), r0d_ptr_a)
440  call mpas_pool_get_array(pool_b, trim(poolitr_b % memberName), r0d_ptr_b)
441  r0d_ptr_b = r0d_ptr_a
442  else if (poolitr_b % nDims == 1) then
443  call mpas_pool_get_array(pool_a, trim(poolitr_a % memberName), r1d_ptr_a)
444  call mpas_pool_get_array(pool_b, trim(poolitr_b % memberName), r1d_ptr_b)
445  r1d_ptr_b = r1d_ptr_a
446  write(message,*) 'Copy sub2all field MIN/MAX: ',trim(poolitr_b % memberName), &
447  minval(r1d_ptr_a),maxval(r1d_ptr_a)
448  call fckit_log%debug(message)
449  else if (poolitr_b % nDims == 2) then
450  call mpas_pool_get_array(pool_a, trim(poolitr_a % memberName), r2d_ptr_a)
451  call mpas_pool_get_array(pool_b, trim(poolitr_b % memberName), r2d_ptr_b)
452  r2d_ptr_b = r2d_ptr_a
453  write(message,*) 'Copy sub2all field MIN/MAX: ',trim(poolitr_b % memberName), &
454  minval(r2d_ptr_a),maxval(r2d_ptr_a)
455  call fckit_log%debug(message)
456  end if
457 
458  else if ( match_scalar(trim(poolitr_b % memberName), trim(poolitr_a % memberName)) ) then
459  write(message,*) 'Copy sub2all field: Looking at SCALARS now',trim(poolitr_a % memberName)
460  call fckit_log%debug(message)
461  call mpas_pool_get_dimension(state, 'index_'//trim(poolitr_a % memberName), index_scalar)
462  if (index_scalar .gt. 0) then
463  call mpas_pool_get_field(pool_a, trim(poolitr_a % memberName), field2d)
464  call mpas_pool_get_field(pool_b, trim(poolitr_b % memberName), field3d)
465  field3d % array(index_scalar,:,:) = field2d % array(:,:)
466  write(message,*) 'Copy sub2all field MIN/MAX: ',trim(poolitr_a % memberName), &
467  minval(field2d % array), maxval(field2d % array)
468  call fckit_log%debug(message)
469  else
470  write(message,*) 'WARNING in Copy sub2all field; ',trim(poolitr_a % memberName), &
471  'not available from MPAS'
472  call fckit_log%debug(message)
473  end if
474 ! end if
475  end if
476  end do
477  end if
478  end if
479  end do
480 
481  end subroutine da_copy_sub2all_fields
482 
483 
484  !***********************************************************************
485  !
486  ! subroutine da_template_pool
487  !
488  !> \brief Subset a pool from fields described in geom
489  !> \details
490  !> Given an mpas_geom object, create templatePool containing a subset
491  !> of the fields in the geom % domain's allFields pool and the
492  !> geom % templated_fields
493  !
494  !-----------------------------------------------------------------------
495  subroutine da_template_pool(geom, templatePool, nf, fieldnames)
496 
497  implicit none
498 
499  ! Arguments
500  type (mpas_geom), intent(in) :: geom
501  type (mpas_pool_type), pointer, intent(out) :: templatepool
502  integer, intent(in) :: nf
503  character (len=*), intent(in) :: fieldnames(nf)
504 
505  ! Local variables
506  type (mpas_pool_type), pointer :: allfields, state
507  character(len=MAXVARLEN) :: fieldname, template
508  !type (field2DReal), pointer :: qField
509  !type (field3DReal), pointer :: scalars
510 
511  integer, pointer :: index_scalar, dim0d
512  integer :: ii
513  integer, parameter :: ndims=10
514  character(len=ShortStrKIND) :: dimnames(ndims)
515 
516  call mpas_pool_create_pool(templatepool)
517 
518  dimnames( 1) = 'nCellsSolve'
519  dimnames( 2) = 'nEdgesSolve'
520  dimnames( 3) = 'nVerticesSolve'
521  dimnames( 4) = 'nVertLevels'
522  dimnames( 5) = 'nVertLevelsP1'
523  dimnames( 6) = 'nSoilLevels'
524  dimnames( 7) = 'nCells'
525  dimnames( 8) = 'nEdges'
526  dimnames( 9) = 'nVertices'
527  dimnames(10) = 'vertexDegree'
528 
529  do ii = 1, ndims
530  call mpas_pool_get_dimension(geom % domain % blocklist % dimensions, trim(dimnames(ii)), dim0d)
531  call mpas_pool_add_dimension(templatepool, trim(dimnames(ii)), dim0d)
532  end do
533 
534  call mpas_pool_get_subpool(geom % domain % blocklist % structs, 'state',state)
535  allfields => geom % domain % blocklist % allFields
536 
537  do ii=1, nf
538  fieldname = fieldnames(ii)
539 
540  if (field_is_scalar(fieldname)) then
541  ! Check if this scalar is activated
542  call mpas_pool_get_dimension(state, 'index_'//trim(fieldname), index_scalar)
543  if (index_scalar .gt. 0) then
544  call mpas_pool_template_field('theta', allfields, fieldname, templatepool)
545  else
546  write(message,*)'--> da_template_pool: ',trim(fieldname), &
547  ' not available in MPAS domain'
548  call abor1_ftn(message)
549  end if
550  else
551  if (pool_has_field(allfields, fieldname)) then
552  template = fieldname
553  else if (geom % is_templated(fieldname)) then
554  template = geom%template(fieldname)
555  if (.not.pool_has_field(allfields, template)) then
556  write(message,*)'--> da_template_pool: ',trim(template), &
557  ' not available in MPAS domain'
558  call abor1_ftn(message)
559  end if
560  else
561  write(message,*)'--> da_template_pool: ',trim(fieldname), &
562  ' not available in MPAS domain or geom % templated_fields'
563  call abor1_ftn(message)
564  end if
565 
566  call mpas_pool_template_field(template, allfields, fieldname, templatepool)
567 
568  end if
569  end do
570 
571  call da_constant(templatepool, mpas_jedi_zero_kr)
572 
573  end subroutine da_template_pool
574 
575 
576  !***********************************************************************
577  !
578  ! subroutine mpas_pool_template_field
579  !
580  !> \brief Add a field to dstPool that is templated on a srcPool field
581  !> \details
582  !> Duplicate srcFieldName from srcPool with the name dstFieldName
583  !> in dstPool
584  !
585  !-----------------------------------------------------------------------
586  subroutine mpas_pool_template_field(srcFieldName, srcPool, dstFieldName, dstPool)
587 
588  ! Arguments
589  type (mpas_pool_type), pointer, intent(in) :: srcPool
590  character (len=*), intent(in) :: srcFieldName, dstFieldName
591  type (mpas_pool_type), pointer, intent(inout) :: dstPool
592 
593  ! Local variables
594  type (mpas_pool_iterator_type) :: poolItr
595  type (field0DReal), pointer :: field0d_src, field0d_dst
596  type (field1DReal), pointer :: field1d_src, field1d_dst
597  type (field2DReal), pointer :: field2d_src, field2d_dst
598  type (field3DReal), pointer :: field3d_src, field3d_dst, field3d
599  type (field0DInteger), pointer :: ifield0d_src, ifield0d_dst
600  type (field1DInteger), pointer :: ifield1d_src, ifield1d_dst
601  type (field2DInteger), pointer :: ifield2d_src, ifield2d_dst
602  type (field3DInteger), pointer :: ifield3d_src, ifield3d_dst
603 
604  !
605  ! Iterate over srcPool, add the one that matches srcFieldName into dstPool
606  !
607  call mpas_pool_begin_iteration(srcpool)
608 
609  do while ( mpas_pool_get_next_member(srcpool, poolitr) )
610 
611  ! Only handle fields, ignore dimensions, namelist options, and other pools
612  if ( poolitr % memberType == mpas_pool_field .and. &
613  trim(srcfieldname) == trim(poolitr % memberName) ) then
614 
615  ! Handle real and integer, ignore logical and char
616  if (poolitr % dataType == mpas_pool_real) then
617 
618  ! Use correctly dimensioned field pointers
619  if (poolitr % nDims == 0) then
620  call mpas_pool_get_field(srcpool, trim(srcfieldname), field0d_src)
621  call mpas_duplicate_field(field0d_src, field0d_dst)
622  field0d_dst % fieldName = dstfieldname
623  call mpas_pool_add_field(dstpool, trim(dstfieldname), field0d_dst)
624  else if (poolitr % nDims == 1) then
625  call mpas_pool_get_field(srcpool, trim(srcfieldname), field1d_src)
626  call mpas_duplicate_field(field1d_src, field1d_dst)
627  field1d_dst % fieldName = dstfieldname
628  call mpas_pool_add_field(dstpool, trim(dstfieldname), field1d_dst)
629  else if (poolitr % nDims == 2) then
630  call mpas_pool_get_field(srcpool, trim(srcfieldname), field2d_src)
631  call mpas_duplicate_field(field2d_src, field2d_dst)
632  field2d_dst % fieldName = dstfieldname
633  call mpas_pool_add_field(dstpool, trim(dstfieldname), field2d_dst)
634  else if (poolitr % nDims == 3) then
635  call mpas_pool_get_field(srcpool, trim(srcfieldname), field3d_src)
636  call mpas_duplicate_field(field3d_src, field3d_dst)
637  field3d_dst % fieldName = dstfieldname
638  call mpas_pool_add_field(dstpool, trim(dstfieldname), field3d_dst)
639  end if
640 
641  else if (poolitr % dataType == mpas_pool_integer) then
642 
643  ! Use correctly dimensioned field pointers
644  if (poolitr % nDims == 0) then
645  call mpas_pool_get_field(srcpool, trim(srcfieldname), ifield0d_src)
646  call mpas_duplicate_field(ifield0d_src, ifield0d_dst)
647  ifield0d_dst % fieldName = dstfieldname
648  call mpas_pool_add_field(dstpool, trim(dstfieldname), ifield0d_dst)
649  else if (poolitr % nDims == 1) then
650  call mpas_pool_get_field(srcpool, trim(srcfieldname), ifield1d_src)
651  call mpas_duplicate_field(ifield1d_src, ifield1d_dst)
652  ifield1d_dst % fieldName = dstfieldname
653  call mpas_pool_add_field(dstpool, trim(dstfieldname), ifield1d_dst)
654  else if (poolitr % nDims == 2) then
655  call mpas_pool_get_field(srcpool, trim(srcfieldname), ifield2d_src)
656  call mpas_duplicate_field(ifield2d_src, ifield2d_dst)
657  ifield2d_dst % fieldName = dstfieldname
658  call mpas_pool_add_field(dstpool, trim(dstfieldname), ifield2d_dst)
659  else if (poolitr % nDims == 3) then
660  call mpas_pool_get_field(srcpool, trim(srcfieldname), ifield3d_src)
661  call mpas_duplicate_field(ifield3d_src, ifield3d_dst)
662  ifield3d_dst % fieldName = dstfieldname
663  call mpas_pool_add_field(dstpool, trim(dstfieldname), ifield3d_dst)
664  end if
665  end if
666 
667  return
668 
669  end if
670  end do
671 
672  write(message,*)'--> mpas_pool_template_field: ',trim(srcfieldname), &
673  ' not available in srcPool'
674  call abor1_ftn(message)
675 
676  end subroutine mpas_pool_template_field
677 
678 
679  !***********************************************************************
680  !
681  ! function da_common_vars
682  !
683  !> \author Gael Descombes
684  !> \date 26 December 2017
685  !> \details
686  !> Count the number of fields in a Pool related to a list of fields
687  !
688  !-----------------------------------------------------------------------
689  function da_common_vars(pool_a, fieldname) result(nsize0)
690 
691  implicit none
692  type (mpas_pool_type), pointer :: pool_a, pool_b
693  type (mpas_pool_iterator_type) :: poolitr
694  character (len=*) :: fieldname(:)
695  integer :: ii, jj, nsize, nsize0
696 
697  nsize0 = 0
698  nsize = size(fieldname)
699  call mpas_pool_begin_iteration(pool_a)
700 
701  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
702  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
703  ! so we select only those members of the pool that are fields
704  if (poolitr % memberType == mpas_pool_field) then
705  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
706  if (poolitr % dataType == mpas_pool_real) then
707  do ii=1, nsize
708  if ( trim(fieldname(ii)).eq.(trim(poolitr % memberName)) ) then
709  nsize0 = nsize0 + 1
710  else if (match_scalar(trim(poolitr % memberName), trim(fieldname(ii)))) then
711  nsize0 = nsize0 + 1
712  end if
713  end do
714  end if
715  end if
716  end do
717 
718  end function da_common_vars
719 
720 
721 
722  !***********************************************************************
723  !
724  ! subroutine da_random
725  !
726  !> \brief Performs random for pool A
727  !> \author Gael Descombes
728  !> \date January 2018
729  !> \details
730  !
731  !-----------------------------------------------------------------------
732  subroutine da_random(pool_a, fld_select)
733 
734  implicit none
735 
736  type (mpas_pool_type), pointer, intent(inout) :: pool_a
737  character (len=*), optional, intent(in) :: fld_select(:)
738 
739  integer, parameter :: rseed = 7
740  type (mpas_pool_iterator_type) :: poolitr
741  real (kind=kind_real), pointer :: r0d_ptr_a
742  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a
743  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a
744  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a
745 
746  !
747  ! Iterate over all fields in pool_b, adding them to fields of the same
748  ! name in pool_a
749  !
750  call mpas_pool_begin_iteration(pool_a)
751 
752  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
753 
754  if (present(fld_select)) then
755  if (ufo_vars_getindex(fld_select,trim(poolitr % memberName)) < 0) cycle
756  end if
757 
758  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
759  ! so we select only those members of the pool that are fields
760  if (poolitr % memberType == mpas_pool_field) then
761 
762  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
763  if (poolitr % dataType == mpas_pool_real) then
764 
765  ! Depending on the dimensionality of the field, we need to set pointers of
766  ! the correct type
767  if (poolitr % nDims == 0) then
768  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
769  !call normal_distribution(r0d_ptr_a, MPAS_JEDI_ZERO_kr, MPAS_JEDI_ONE_kr, rseed)
770  else if (poolitr % nDims == 1) then
771  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
772  call normal_distribution(r1d_ptr_a, mpas_jedi_zero_kr, mpas_jedi_one_kr, rseed)
773  else if (poolitr % nDims == 2) then
774  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
775  call normal_distribution(r2d_ptr_a, mpas_jedi_zero_kr, mpas_jedi_one_kr, rseed)
776  else if (poolitr % nDims == 3) then
777  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
778  call normal_distribution(r3d_ptr_a, mpas_jedi_zero_kr, mpas_jedi_one_kr, rseed)
779  end if
780 
781  end if
782  end if
783  end do
784 
785  end subroutine da_random
786 
787  !-----------------------------------------------------------------------
788  ! subroutine da_operator
789  !
790  !> \brief Performs A = A 'kind_op' B for pools A and B
791  !> \author Michael Duda
792  !> \date 20 December 2017
793  !> \details
794  !> Given two pools, A and B, where the fields in B are a subset of
795  !> the fields in A, this routine adds the fields in B to fields in A
796  !> with the same name. When A and B contain identical fields, this
797  !> is equivalent to A = A 'kind_op' B.
798  !> \modified by Gael DESCOMBES to apply diffferent operator
799  !
800  !-----------------------------------------------------------------------
801  subroutine da_operator(kind_op, pool_a, pool_b, pool_c, fld_select)
802 
803  implicit none
804 
805  type (mpas_pool_type), pointer :: pool_a, pool_b
806  type (mpas_pool_type), pointer, optional :: pool_c
807  character (len=*) :: kind_op
808  character (len=*), optional :: fld_select(:)
809 
810  type (mpas_pool_iterator_type) :: poolitr
811  real (kind=kind_real), pointer :: r0d_ptr_a, r0d_ptr_b, r0d_ptr_c
812  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a, r1d_ptr_b, r1d_ptr_c
813  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a, r2d_ptr_b, r2d_ptr_c
814  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a, r3d_ptr_b, r3d_ptr_c
815  !
816  ! Iterate over all fields in pool_b, adding them to fields of the same
817  ! name in pool_a
818  !
819  call mpas_pool_begin_iteration(pool_b)
820 
821  write(message,*) 'Operator ',trim(kind_op)
822  call fckit_log%debug(message)
823 
824  do while ( mpas_pool_get_next_member(pool_b, poolitr) )
825 
826  if (present(fld_select)) then
827  if (ufo_vars_getindex(fld_select,trim(poolitr % memberName)) < 0) cycle
828  end if
829 
830  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
831  ! so we select only those members of the pool that are fields
832  if (poolitr % memberType == mpas_pool_field) then
833 
834  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
835  if (poolitr % dataType == mpas_pool_real) then
836 
837  ! Depending on the dimensionality of the field, we need to set pointers of
838  ! the correct type
839  if (poolitr % nDims == 0) then
840  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
841  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r0d_ptr_b)
842  if (present(pool_c)) then
843  call mpas_pool_get_array(pool_c, trim(poolitr % memberName), r0d_ptr_c)
844  r0d_ptr_a = mpas_jedi_zero_kr
845  end if
846  if ( trim(kind_op).eq.'add' ) then
847  r0d_ptr_a = r0d_ptr_a + r0d_ptr_b
848  if (present(pool_c)) then
849  r0d_ptr_a = r0d_ptr_b + r0d_ptr_c
850  else
851  r0d_ptr_a = r0d_ptr_a + r0d_ptr_b
852  end if
853  else if ( trim(kind_op).eq.'schur' ) then
854  if (present(pool_c)) then
855  r0d_ptr_a = r0d_ptr_b * r0d_ptr_c
856  else
857  r0d_ptr_a = r0d_ptr_a * r0d_ptr_b
858  end if
859  else if ( trim(kind_op).eq.'sub' ) then
860  if (present(pool_c)) then
861  r0d_ptr_a = r0d_ptr_b - r0d_ptr_c
862  else
863  r0d_ptr_a = r0d_ptr_a - r0d_ptr_b
864  end if
865  end if
866 
867  else if (poolitr % nDims == 1) then
868  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
869  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r1d_ptr_b)
870  if (present(pool_c)) then
871  call mpas_pool_get_array(pool_c, trim(poolitr % memberName), r1d_ptr_c)
872  r1d_ptr_a = mpas_jedi_zero_kr
873  end if
874  if ( trim(kind_op).eq.'add' ) then
875  if (present(pool_c)) then
876  r1d_ptr_a = r1d_ptr_b + r1d_ptr_c
877  else
878  r1d_ptr_a = r1d_ptr_a + r1d_ptr_b
879  end if
880  else if ( trim(kind_op).eq.'schur' ) then
881  if (present(pool_c)) then
882  r1d_ptr_a = r1d_ptr_b * r1d_ptr_c
883  else
884  r1d_ptr_a = r1d_ptr_a * r1d_ptr_b
885  end if
886  else if ( trim(kind_op).eq.'sub' ) then
887  if (present(pool_c)) then
888  r1d_ptr_a = r1d_ptr_b - r1d_ptr_c
889  else
890  r1d_ptr_a = r1d_ptr_a - r1d_ptr_b
891  end if
892  end if
893 
894  else if (poolitr % nDims == 2) then
895  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
896  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r2d_ptr_b)
897  if (present(pool_c)) then
898  call mpas_pool_get_array(pool_c, trim(poolitr % memberName), r2d_ptr_c)
899  r2d_ptr_a = mpas_jedi_zero_kr
900  end if
901  if ( trim(kind_op).eq.'add' ) then
902  write(message,*) 'Operator_a add MIN/MAX: ',minval(r2d_ptr_a),maxval(r2d_ptr_a)
903  call fckit_log%debug(message)
904  write(message,*) 'Operator_b add MIN/MAX: ',minval(r2d_ptr_b),maxval(r2d_ptr_b)
905  call fckit_log%debug(message)
906  if (present(pool_c)) then
907  r2d_ptr_a = r2d_ptr_b + r2d_ptr_c
908  else
909  r2d_ptr_a = r2d_ptr_a + r2d_ptr_b
910  end if
911  write(message,*) 'Operator2 add MIN/MAX: ',minval(r2d_ptr_a),maxval(r2d_ptr_a)
912  call fckit_log%debug(message)
913  else if ( trim(kind_op).eq.'schur' ) then
914  if (present(pool_c)) then
915  r2d_ptr_a = r2d_ptr_b * r2d_ptr_c
916  else
917  r2d_ptr_a = r2d_ptr_a * r2d_ptr_b
918  end if
919  else if ( trim(kind_op).eq.'sub' ) then
920  if (present(pool_c)) then
921  r2d_ptr_a = r2d_ptr_b - r2d_ptr_c
922  else
923  r2d_ptr_a = r2d_ptr_a - r2d_ptr_b
924  end if
925  end if
926 
927  else if (poolitr % nDims == 3) then
928  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
929  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r3d_ptr_b)
930  if (present(pool_c)) then
931  call mpas_pool_get_array(pool_c, trim(poolitr % memberName), r3d_ptr_c)
932  r3d_ptr_a = mpas_jedi_zero_kr
933  end if
934  if ( trim(kind_op).eq.'add' ) then
935  if (present(pool_c)) then
936  r3d_ptr_a = r3d_ptr_b + r3d_ptr_c
937  else
938  r3d_ptr_a = r3d_ptr_a + r3d_ptr_b
939  end if
940  else if ( trim(kind_op).eq.'schur' ) then
941  if (present(pool_c)) then
942  r3d_ptr_a = r3d_ptr_b * r3d_ptr_c
943  else
944  r3d_ptr_a = r3d_ptr_a * r3d_ptr_b
945  end if
946  else if ( trim(kind_op).eq.'sub' ) then
947  if (present(pool_c)) then
948  r3d_ptr_a = r3d_ptr_b - r3d_ptr_c
949  else
950  r3d_ptr_a = r3d_ptr_a - r3d_ptr_b
951  end if
952  end if
953  end if
954 
955  end if
956  end if
957  end do
958 
959  end subroutine da_operator
960 
961  !***********************************************************************
962  !
963  ! subroutine da_self_mult
964  !
965  !> \brief Performs A = A * zz for pool A, zz a real number
966  !> \author Gael Descombes
967  !> \date 22 December 2017
968  !> \details
969  !
970  !-----------------------------------------------------------------------
971  subroutine da_self_mult(pool_a, zz)
972 
973  implicit none
974 
975  type (mpas_pool_type), pointer :: pool_a
976  real (kind=kind_real) :: zz
977 
978  type (mpas_pool_iterator_type) :: poolitr
979  real (kind=kind_real), pointer :: r0d_ptr_a
980  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a
981  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a
982  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a
983 
984  !
985  ! Iterate over all fields in pool_b, adding them to fields of the same
986  ! name in pool_a
987  !
988  call mpas_pool_begin_iteration(pool_a)
989 
990  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
991 
992  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
993  ! so we select only those members of the pool that are fields
994  if (poolitr % memberType == mpas_pool_field) then
995 
996  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
997  if (poolitr % dataType == mpas_pool_real) then
998 
999  ! Depending on the dimensionality of the field, we need to set pointers of
1000  ! the correct type
1001  if (poolitr % nDims == 0) then
1002  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
1003  r0d_ptr_a = r0d_ptr_a * zz
1004  else if (poolitr % nDims == 1) then
1005  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
1006  r1d_ptr_a = r1d_ptr_a * zz
1007  else if (poolitr % nDims == 2) then
1008  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
1009  r2d_ptr_a = r2d_ptr_a * zz
1010  else if (poolitr % nDims == 3) then
1011  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
1012  r3d_ptr_a = r3d_ptr_a * zz
1013  end if
1014 
1015  end if
1016  end if
1017  end do
1018 
1019  end subroutine da_self_mult
1020 
1021 
1022  !***********************************************************************
1023  !
1024  ! subroutine da_constant
1025  !
1026  !> \brief Performs A = constant. for pool A
1027  !> \author Gael Descombes
1028  !> \date 22 December 2017
1029  !> \details
1030  !
1031  !-----------------------------------------------------------------------
1032  subroutine da_constant(pool_a, realvalue, fld_select)
1033 
1034  implicit none
1035 
1036  type (mpas_pool_type), pointer, intent(inout) :: pool_a
1037  real (kind=kind_real), intent(in) :: realvalue
1038  character (len=*), optional, intent(in) :: fld_select(:)
1039 
1040  type (mpas_pool_iterator_type) :: poolitr
1041  real (kind=kind_real), pointer :: r0d_ptr_a
1042  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a
1043  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a
1044  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a
1045 
1046  call mpas_pool_begin_iteration(pool_a)
1047 
1048  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
1049 
1050  if (present(fld_select)) then
1051  if (ufo_vars_getindex(fld_select,trim(poolitr % memberName)) < 0) cycle
1052  end if
1053 
1054  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
1055  ! so we select only those members of the pool that are fields
1056  if (poolitr % memberType == mpas_pool_field) then
1057 
1058  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
1059  if (poolitr % dataType == mpas_pool_real) then
1060 
1061  ! Depending on the dimensionality of the field, we need to set pointers of
1062  ! the correct type
1063  if (poolitr % nDims == 0) then
1064  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
1065  r0d_ptr_a = realvalue
1066  else if (poolitr % nDims == 1) then
1067  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
1068  r1d_ptr_a = realvalue
1069  else if (poolitr % nDims == 2) then
1070  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
1071  r2d_ptr_a = realvalue
1072  else if (poolitr % nDims == 3) then
1073  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
1074  r3d_ptr_a = realvalue
1075  end if
1076 
1077  end if
1078  end if
1079  end do
1080 
1081  end subroutine da_constant
1082 
1083  !***********************************************************************
1084  !
1085  ! subroutine da_posdef
1086  !
1087  !> \brief Performs A = max(0.,A) for pool A
1088  !> \author JJ Guerrette
1089  !> \date 12 July 2019
1090  !> \details
1091  !
1092  !-----------------------------------------------------------------------
1093  subroutine da_posdef(pool_a, fld_select)
1094 
1095  implicit none
1096 
1097  type (mpas_pool_type), pointer, intent(inout) :: pool_a
1098  character (len=*), optional, intent(in) :: fld_select(:)
1099 
1100  type (mpas_pool_iterator_type) :: poolitr
1101  real (kind=kind_real), pointer :: r0d_ptr_a
1102  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a
1103  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a
1104  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a
1105 
1106  !
1107  ! Iterate over all fields in pool_b, adding them to fields of the same
1108  ! name in pool_a
1109  !
1110  call mpas_pool_begin_iteration(pool_a)
1111 
1112  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
1113 
1114  if (present(fld_select)) then
1115  if (ufo_vars_getindex(fld_select,trim(poolitr % memberName)) < 0) cycle
1116  end if
1117 
1118  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
1119  ! so we select only those members of the pool that are fields
1120  if (poolitr % memberType == mpas_pool_field) then
1121 
1122  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
1123  if (poolitr % dataType == mpas_pool_real) then
1124 
1125  ! Depending on the dimensionality of the field, we need to set pointers of
1126  ! the correct type
1127  if (poolitr % nDims == 0) then
1128  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
1129  r0d_ptr_a = max(mpas_jedi_zero_kr, r0d_ptr_a)
1130  else if (poolitr % nDims == 1) then
1131  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
1132  r1d_ptr_a = max(mpas_jedi_zero_kr, r1d_ptr_a)
1133  else if (poolitr % nDims == 2) then
1134  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
1135  r2d_ptr_a = max(mpas_jedi_zero_kr, r2d_ptr_a)
1136  else if (poolitr % nDims == 3) then
1137  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
1138  r3d_ptr_a = max(mpas_jedi_zero_kr, r3d_ptr_a)
1139  end if
1140 
1141  end if
1142  end if
1143  end do
1144 
1145  end subroutine da_posdef
1146 
1147  !***********************************************************************
1148  !
1149  ! subroutine da_setval
1150  !
1151  !> \brief Performs A = Val_R. for pool A
1152  !> \author Gael Descombes
1153  !> \date 22 December 2017
1154  !> \details
1155  !
1156  !-----------------------------------------------------------------------
1157  subroutine da_setval(pool_a,zz)
1158 
1159  implicit none
1160 
1161  type (mpas_pool_type), pointer :: pool_a
1162  real (kind=kind_real) :: zz
1163 
1164  type (mpas_pool_iterator_type) :: poolitr
1165  real (kind=kind_real), pointer :: r0d_ptr_a
1166  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a
1167  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a
1168  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a
1169 
1170  !
1171  ! Iterate over all fields in pool_b, adding them to fields of the same
1172  ! name in pool_a
1173  !
1174  call mpas_pool_begin_iteration(pool_a)
1175 
1176  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
1177 
1178  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
1179  ! so we select only those members of the pool that are fields
1180  if (poolitr % memberType == mpas_pool_field) then
1181 
1182  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
1183  if (poolitr % dataType == mpas_pool_real) then
1184 
1185  ! Depending on the dimensionality of the field, we need to set pointers of
1186  ! the correct type
1187  if (poolitr % nDims == 0) then
1188  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
1189  r0d_ptr_a = zz
1190  else if (poolitr % nDims == 1) then
1191  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
1192  r1d_ptr_a = zz
1193  else if (poolitr % nDims == 2) then
1194  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
1195  r2d_ptr_a = zz
1196  else if (poolitr % nDims == 3) then
1197  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
1198  r3d_ptr_a = zz
1199  end if
1200 
1201  end if
1202  end if
1203  end do
1204 
1205  end subroutine da_setval
1206 
1207 
1208 
1209 
1210 
1211  !***********************************************************************
1212  !
1213  ! subroutine da_axpy
1214  !
1215  !> \brief Performs A = A + B * zz for pools A and B
1216  !> \author Gael Descombes
1217  !> \date 20 December 2017
1218  !> \details
1219  !> Given two pools, A and B, where the fields in B are a subset of
1220  !> the fields in A, this routine adds the fields in B to fields in A
1221  !> with the same name. When A and B contain identical fields, this
1222  !> is equivalent to A = A + B.
1223  !
1224  !-----------------------------------------------------------------------
1225  subroutine da_axpy(pool_a, pool_b, zz, fld_select)
1226 
1227  implicit none
1228  type (mpas_pool_type), pointer, intent(inout) :: pool_a
1229  type (mpas_pool_type), pointer, intent(in) :: pool_b
1230  real (kind=kind_real), intent(in) :: zz
1231  character (len=*), optional, intent(in) :: fld_select(:)
1232 
1233 
1234 
1235  type (mpas_pool_iterator_type) :: poolitr
1236  real (kind=kind_real), pointer :: r0d_ptr_a, r0d_ptr_b
1237  real (kind=kind_real), dimension(:), pointer :: r1d_ptr_a, r1d_ptr_b
1238  real (kind=kind_real), dimension(:,:), pointer :: r2d_ptr_a, r2d_ptr_b
1239  real (kind=kind_real), dimension(:,:,:), pointer :: r3d_ptr_a, r3d_ptr_b
1240 
1241  !
1242  ! Iterate over all fields in pool_b, adding them to fields of the same
1243  ! name in pool_a
1244  !
1245  call mpas_pool_begin_iteration(pool_b)
1246 
1247  do while ( mpas_pool_get_next_member(pool_b, poolitr) )
1248 
1249  if (present(fld_select)) then
1250  if (ufo_vars_getindex(fld_select,trim(poolitr % memberName)) < 0) cycle
1251  end if
1252 
1253  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
1254  ! so we select only those members of the pool that are fields
1255  if (poolitr % memberType == mpas_pool_field) then
1256 
1257  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
1258  if (poolitr % dataType == mpas_pool_real) then
1259 
1260  ! Depending on the dimensionality of the field, we need to set pointers of
1261  ! the correct type
1262  if (poolitr % nDims == 0) then
1263  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r0d_ptr_a)
1264  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r0d_ptr_b)
1265  r0d_ptr_a = r0d_ptr_a + r0d_ptr_b * zz
1266  else if (poolitr % nDims == 1) then
1267  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r1d_ptr_a)
1268  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r1d_ptr_b)
1269  r1d_ptr_a = r1d_ptr_a + r1d_ptr_b * zz
1270  else if (poolitr % nDims == 2) then
1271  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r2d_ptr_a)
1272  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r2d_ptr_b)
1273  r2d_ptr_a = r2d_ptr_a + r2d_ptr_b * zz
1274  else if (poolitr % nDims == 3) then
1275  call mpas_pool_get_array(pool_a, trim(poolitr % memberName), r3d_ptr_a)
1276  call mpas_pool_get_array(pool_b, trim(poolitr % memberName), r3d_ptr_b)
1277  r3d_ptr_a = r3d_ptr_a + r3d_ptr_b * zz
1278  end if
1279 
1280  end if
1281  end if
1282  end do
1283 
1284  end subroutine da_axpy
1285 
1286 
1287  !***********************************************************************
1288  !
1289  ! subroutine da_gpnorm
1290  !
1291  !> \brief Performs basic statistics min/max/norm given a pool
1292  !> \author Gael Descombes
1293  !> \date February 2018
1294  !> \details
1295  !> Given a pool of fields, return min/max/norm array
1296  !
1297  !-----------------------------------------------------------------------
1298 
1299  subroutine da_gpnorm(pool_a, dminfo, nf, pstat, fld_select)
1300 
1301  implicit none
1302  type (mpas_pool_type), pointer, intent(in) :: pool_a
1303  type (dm_info), pointer, intent(in) :: dminfo
1304  integer, intent(in) :: nf
1305  character (len=*), intent(in) :: fld_select(nf)
1306  real(kind=kind_real), intent(out) :: pstat(3, nf)
1307 
1308  type (mpas_pool_iterator_type) :: poolitr
1309  type (field1dreal), pointer :: field1d
1310  type (field2dreal), pointer :: field2d
1311  type (field3dreal), pointer :: field3d
1312  real(kind=kind_real) :: globalsum, globalmin, globalmax, dimtot, dimtot_global, prodtot
1313 
1314  integer :: jj, ndims
1315  integer :: dim1, dim2, dim3
1316  integer, allocatable :: dimsizes(:)
1317 
1318  pstat = mpas_jedi_zero_kr
1319 
1320  !
1321  ! Iterate over all fields in pool_a
1322  ! name in pool_a
1323  !
1324  call mpas_pool_begin_iteration(pool_a)
1325 
1326  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
1327  jj = ufo_vars_getindex(fld_select,trim(poolitr % memberName))
1328  if ( jj < 0 .or. jj > nf ) cycle
1329 
1330  ! Pools may in general contain dimensions, namelist options, fields, or other pools,
1331  ! so we select only those members of the pool that are fields
1332  if (poolitr % memberType == mpas_pool_field) then
1333 
1334  ! Fields can be integer, logical, or real. Here, we operate only on real-valued fields
1335  if (poolitr % dataType == mpas_pool_real) then
1336 
1337  ! Depending on the dimensionality of the field, we need to set pointers of
1338  ! the correct type
1339  ndims = poolitr % nDims
1340  dimsizes = getsolvedimsizes(pool_a, poolitr%memberName)
1341  if (ndims == 1) then
1342  dim1 = dimsizes(1)
1343  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field1d)
1344  dimtot = real(dim1,kind_real)
1345  prodtot = sum(field1d % array(1:dim1)**2 )
1346  call mpas_dmpar_sum_real(dminfo, dimtot, dimtot_global)
1347  call mpas_dmpar_sum_real(dminfo, prodtot, globalsum)
1348  call mpas_dmpar_min_real(dminfo, minval(field1d % array(1:dim1)), globalmin)
1349  call mpas_dmpar_max_real(dminfo, maxval(field1d % array(1:dim1)), globalmax)
1350  pstat(1,jj) = globalmin
1351  pstat(2,jj) = globalmax
1352  pstat(3,jj) = sqrt( globalsum / dimtot_global )
1353  else if (ndims == 2) then
1354  dim1 = dimsizes(1)
1355  dim2 = dimsizes(2)
1356  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field2d)
1357  dimtot = real(dim1*dim2,kind_real)
1358  prodtot = sum(field2d % array(1:dim1,1:dim2)**2 )
1359  call mpas_dmpar_sum_real(dminfo, dimtot, dimtot_global)
1360  call mpas_dmpar_sum_real(dminfo, prodtot, globalsum)
1361  call mpas_dmpar_min_real(dminfo, minval(field2d % array(1:dim1,1:dim2)), globalmin)
1362  call mpas_dmpar_max_real(dminfo, maxval(field2d % array(1:dim1,1:dim2)), globalmax)
1363  pstat(1,jj) = globalmin
1364  pstat(2,jj) = globalmax
1365  pstat(3,jj) = sqrt( globalsum / dimtot_global )
1366  else if (ndims == 3) then
1367  dim1 = dimsizes(1)
1368  dim2 = dimsizes(2)
1369  dim3 = dimsizes(3)
1370  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field3d)
1371  dimtot = real(dim1*dim2*dim3,kind_real)
1372  prodtot = sum(field3d % array(1:dim1,1:dim2,1:dim3)**2 )
1373  call mpas_dmpar_sum_real(dminfo, dimtot, dimtot_global)
1374  call mpas_dmpar_sum_real(dminfo, prodtot, globalsum)
1375  call mpas_dmpar_min_real(dminfo, minval(field3d % array(1:dim1,1:dim2,1:dim3)), globalmin)
1376  call mpas_dmpar_max_real(dminfo, maxval(field3d % array(1:dim1,1:dim2,1:dim3)), globalmax)
1377  pstat(1,jj) = globalmin
1378  pstat(2,jj) = globalmax
1379  pstat(3,jj) = sqrt( globalsum / dimtot_global )
1380  end if
1381  deallocate(dimsizes)
1382  end if
1383  end if
1384  end do
1385 
1386  end subroutine da_gpnorm
1387 
1388 
1389  !***********************************************************************
1390  !
1391  ! subroutine da_fldrms
1392  !
1393  !> \brief Performs basic statistics min/max/norm given a pool
1394  !> \author Gael Descombes
1395  !> \date February 2018
1396  !> \details
1397  !> Given a pool of fields, return min/max/norm array
1398  !
1399  !-----------------------------------------------------------------------
1400 
1401  subroutine da_fldrms(pool_a, dminfo, fldrms, fld_select)
1402 
1403  implicit none
1404  type (mpas_pool_type), pointer, intent(in) :: pool_a
1405  type (dm_info), pointer, intent(in) :: dminfo
1406  real(kind=kind_real), intent(out) :: fldrms
1407  character (len=*), optional, intent(in) :: fld_select(:)
1408 
1409  type (mpas_pool_iterator_type) :: poolitr
1410  type (field1dreal), pointer :: field1d
1411  type (field2dreal), pointer :: field2d
1412  type (field3dreal), pointer :: field3d
1413  real(kind=kind_real) :: dimtot, dimtot_global, prodtot, prodtot_global
1414 
1415  integer :: ndims
1416  integer :: dim1, dim2, dim3
1417  integer, allocatable :: dimsizes(:)
1418 
1419  prodtot = mpas_jedi_zero_kr
1420  dimtot = mpas_jedi_zero_kr
1421 
1422  !
1423  ! Iterate over all fields in pool_a
1424  ! named in pool_a
1425  !
1426  call mpas_pool_begin_iteration(pool_a)
1427 
1428  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
1429  if (present(fld_select)) then
1430  if (ufo_vars_getindex(fld_select,trim(poolitr % memberName)) < 0) cycle
1431  end if
1432  if (poolitr % dataType == mpas_pool_real) then
1433  if (poolitr % memberType == mpas_pool_field) then
1434  ndims = poolitr % nDims
1435  dimsizes = getsolvedimsizes(pool_a, poolitr%memberName)
1436  if (ndims == 1) then
1437  dim1 = dimsizes(1)
1438  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field1d)
1439  dimtot = dimtot + real(dim1,kind_real)
1440  prodtot = prodtot + sum( field1d % array(1:dim1)**2 )
1441  else if (ndims == 2) then
1442  dim1 = dimsizes(1)
1443  dim2 = dimsizes(2)
1444  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field2d)
1445  dimtot = dimtot + real(dim1*dim2,kind_real)
1446  prodtot = prodtot + sum( field2d % array(1:dim1,1:dim2)**2 )
1447  else if (ndims == 3) then
1448  dim1 = dimsizes(1)
1449  dim2 = dimsizes(2)
1450  dim3 = dimsizes(3)
1451  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field3d)
1452  dimtot = dimtot + real(dim1*dim2*dim3,kind_real)
1453  prodtot = prodtot + sum( field3d % array(1:dim1,1:dim2,1:dim3)**2 )
1454  end if
1455  deallocate(dimsizes)
1456  end if
1457  end if
1458  end do
1459 
1460  call mpas_dmpar_sum_real(dminfo, dimtot, dimtot_global)
1461  call mpas_dmpar_sum_real(dminfo, prodtot, prodtot_global)
1462  fldrms = sqrt(prodtot_global / dimtot_global)
1463 
1464  end subroutine da_fldrms
1465 
1466 
1467  !***********************************************************************
1468  !
1469  ! subroutine da_dot_product
1470  !
1471  !> \brief Performs the dot_product given two pools of fields
1472  !> \author Gael Descombes
1473  !> \date February 2018
1474  !> \details
1475  !> Given two pools of fields, compute the dot_product
1476  !
1477  !-----------------------------------------------------------------------
1478 
1479  subroutine da_dot_product(pool_a, pool_b, dminfo, zprod)
1480 
1481  implicit none
1482  type (mpas_pool_type), pointer, intent(in) :: pool_a, pool_b
1483  type (dm_info), pointer, intent(in) :: dminfo
1484  real(kind=kind_real), intent(out) :: zprod
1485 
1486  type (mpas_pool_iterator_type) :: poolitr
1487  type (field1dreal), pointer :: field1d_a, field1d_b
1488  type (field2dreal), pointer :: field2d_a, field2d_b
1489  type (field3dreal), pointer :: field3d_a, field3d_b
1490  real(kind=kind_real) :: fieldsum_local, zprod_local
1491 
1492  integer :: ndims
1493  integer :: dim1, dim2, dim3
1494  integer, allocatable :: dimsizes(:)
1495 
1496  !
1497  ! Iterate over all fields in pool_a
1498  ! named in pool_a
1499  !
1500  call mpas_pool_begin_iteration(pool_a)
1501 
1502  zprod_local = mpas_jedi_zero_kr
1503 
1504  do while ( mpas_pool_get_next_member(pool_a, poolitr) )
1505  if (poolitr % dataType == mpas_pool_real) then
1506  if (poolitr % memberType == mpas_pool_field) then
1507  ndims = poolitr % nDims
1508  dimsizes = getsolvedimsizes(pool_a, poolitr%memberName)
1509  !TODO: add check that dimSizes are the same between pool_a and pool_b for poolItr
1510  if (ndims == 1) then
1511  dim1 = dimsizes(1)
1512  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field1d_a)
1513  call mpas_pool_get_field(pool_b, trim(poolitr % memberName), field1d_b)
1514  fieldsum_local = sum(field1d_a % array(1:dim1) * field1d_b % array(1:dim1))
1515  zprod_local = zprod_local + fieldsum_local
1516  else if (ndims == 2) then
1517  dim1 = dimsizes(1)
1518  dim2 = dimsizes(2)
1519  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field2d_a)
1520  call mpas_pool_get_field(pool_b, trim(poolitr % memberName), field2d_b)
1521  fieldsum_local = sum(field2d_a % array(1:dim1,1:dim2) &
1522  * field2d_b % array(1:dim1,1:dim2))
1523  zprod_local = zprod_local + fieldsum_local
1524  else if (ndims == 3) then
1525  dim1 = dimsizes(1)
1526  dim2 = dimsizes(2)
1527  dim3 = dimsizes(3)
1528  call mpas_pool_get_field(pool_a, trim(poolitr % memberName), field3d_a)
1529  call mpas_pool_get_field(pool_b, trim(poolitr % memberName), field3d_b)
1530  fieldsum_local = sum(field3d_a % array(1:dim1,1:dim2,1:dim3) &
1531  * field3d_b % array(1:dim1,1:dim2,1:dim3))
1532  zprod_local = zprod_local + fieldsum_local
1533  end if
1534  deallocate(dimsizes)
1535  end if
1536  end if
1537  end do
1538 
1539  call mpas_dmpar_sum_real(dminfo, zprod_local, zprod)
1540 
1541  end subroutine da_dot_product
1542 
1543 
1544  subroutine cvt_oopsmpas_date(inString2,outString2,iconv)
1545 
1546  implicit none
1547 
1548  character (len=*), intent(in) :: instring2
1549  character (len=*), intent(inout) :: outstring2
1550  integer, intent(in) :: iconv
1551  integer :: i, curlen
1552  integer :: year, month, day, hour, minute, second
1553 
1554  character (len=ShortStrKIND) :: timepart
1555  character (len=ShortStrKIND) :: yearformat
1556  logical :: charexpand
1557  character (len=4) :: yyyy
1558  character (len=2) :: mm, dd, h, m, s
1559  character (len=21) :: outstring, instring
1560 
1561  ! 2017-08-08T00:00:00Z OOPS/YAML format
1562  ! 2010-10-24_02.00.00 MPAS format
1563  ! iconv=1: MPAS --> OOPS/YAML
1564  ! iconv=-1: OOPS/YAML --> MPAS
1565 
1566  if (iconv.eq.-1) then
1567  yyyy = instring2(1:4)
1568  mm = instring2(6:7)
1569  dd = instring2(9:10)
1570  h = instring2(12:13)
1571  m = instring2(15:16)
1572  s = instring2(18:19)
1573  else
1574  yyyy = instring2(1:4)
1575  mm = instring2(6:7)
1576  dd = instring2(9:10)
1577  h = instring2(12:13)
1578  m = instring2(15:16)
1579  s = instring2(18:19)
1580  end if
1581 
1582  write(message,*) 'cvt_oopsmpas_date instring: ',trim(yyyy),trim(mm),trim(dd),trim(h),trim(m),trim(s)
1583  call fckit_log%debug(message)
1584  write(message,*) 'cvt_oopsmpas_date input ',trim(instring2)
1585  call fckit_log%debug(message)
1586 
1587  write(outstring,*) ''
1588  instring = trim(outstring2)
1589 
1590  curlen = 0
1591  charexpand = .false.
1592  do i = 1, len_trim(instring)
1593  if (instring(i:i) == '$' ) then
1594  charexpand = .true.
1595  else if (instring(i:i) /= '$') then
1596  write(message,*) 'inString: ',trim(instring(i:i)),charexpand
1597  call fckit_log%debug(message)
1598  if (charexpand) then
1599  select case (instring(i:i))
1600  case ('Y')
1601  outstring = trim(outstring) // trim(yyyy)
1602  case ('M')
1603  outstring = trim(outstring) // trim(mm)
1604  case ('D')
1605  outstring = trim(outstring) // trim(dd)
1606  case ('h')
1607  outstring = trim(outstring) // trim(h)
1608  case ('m')
1609  outstring = trim(outstring) // trim(m)
1610  case ('s')
1611  outstring = trim(outstring) // trim(s)
1612  case default
1613  call mpas_dmpar_global_abort('ERROR: mpas_timekeeping')
1614  end select
1615  curlen = len_trim(outstring)
1616  charexpand = .false.
1617  write(message,*) 'outString: ',trim(outstring)
1618  call fckit_log%debug(message)
1619  else
1620  outstring(curlen+1:curlen+1) = outstring2(i:i)
1621  curlen = curlen+1
1622  end if
1623  end if
1624  end do
1625 
1626  outstring2 = trim(outstring)
1627  write(message,*) 'cvt_oopsmpas_date output ',trim(outstring2)
1628  call fckit_log%debug(message)
1629 
1630  end subroutine cvt_oopsmpas_date
1631 
1632 
1633 
1634 ! ------------------------------------------------------------------------
1635 ! chunk of code from DART
1636 ! ------------------------------------------------------------------------
1637 
1638 subroutine uv_cell_to_edges(domain, u_field, v_field, du, lonCell, latCell, &
1639  & nCells, edgeNormalVectors, nEdgesOnCell, edgesOnCell, nVertLevels)
1640 
1641  ! Project u, v wind increments at cell centers onto the edges.
1642  ! FIXME:
1643  ! we can hard-code R3 here since it comes from the (3d) x/y/z cartesian coordinate.
1644  ! We define nEdgesOnCell in get_grid_dims, and read edgesOnCell in get_grid.
1645  ! We read edgeNormalVectors in get_grid to use this subroutine.
1646  ! Here "U" is the prognostic variable in MPAS, and we update it with the wind
1647  ! increments at cell centers.
1648 
1649  implicit none
1650 
1651  type (domain_type), pointer, intent(inout) :: domain
1652  type (field2dreal), pointer, intent(in) :: u_field ! u wind updated from filter
1653  type (field2dreal), pointer, intent(in) :: v_field ! v wind updated from filter
1654  type (field2dreal), pointer, intent(inout) :: du ! normal velocity increment on the edges
1655  real(kind_real), intent(in) :: loncell(1:ncells), latcell(1:ncells) ! lon, lat at cell centers in radians
1656  real(kind_real), intent(in) :: edgenormalvectors(:,:)
1657  integer, intent(in) :: nedgesoncell(:), edgesoncell(:,:)
1658  integer, intent(in) :: ncells, nvertlevels
1659 
1660  ! Local variables
1661  integer, parameter :: r3 = 3
1662  real(kind_real), dimension(:,:), allocatable :: east, north
1663  integer :: icell, iedge, jedge, k
1664 
1665  ! allocation
1666  allocate(east(r3,ncells))
1667  allocate(north(r3,ncells))
1668 
1669  ! Initialization
1670  du%array(:,:) = mpas_jedi_zero_kr
1671 
1672  ! Compute unit vectors in east and north directions for each cell:
1673  do icell = 1, ncells
1674  east(1,icell) = -sin(loncell(icell))
1675  east(2,icell) = cos(loncell(icell))
1676  east(3,icell) = mpas_jedi_zero_kr
1677  call r3_normalize(east(1,icell), east(2,icell), east(3,icell))
1678  north(1,icell) = -cos(loncell(icell))*sin(latcell(icell))
1679  north(2,icell) = -sin(loncell(icell))*sin(latcell(icell))
1680  north(3,icell) = cos(latcell(icell))
1681  call r3_normalize(north(1,icell), north(2,icell), north(3,icell))
1682  end do
1683 
1684  ! Project analysis increments from the cell centers to the edges
1685 
1686  do icell = 1, ncells
1687  do jedge = 1, nedgesoncell(icell)
1688  iedge = edgesoncell(jedge, icell)
1689  do k = 1, nvertlevels
1690  du%array(k,iedge) = du%array(k,iedge) + mpas_jedi_half_kr * u_field%array(k,icell) &
1691  * (edgenormalvectors(1,iedge) * east(1,icell) &
1692  + edgenormalvectors(2,iedge) * east(2,icell) &
1693  + edgenormalvectors(3,iedge) * east(3,icell)) &
1694  + mpas_jedi_half_kr * v_field%array(k,icell) &
1695  * (edgenormalvectors(1,iedge) * north(1,icell) &
1696  + edgenormalvectors(2,iedge) * north(2,icell) &
1697  + edgenormalvectors(3,iedge) * north(3,icell))
1698  end do
1699  end do
1700  end do
1701 
1702  ! deallocation
1703  deallocate(east)
1704  deallocate(north)
1705 
1706 end subroutine uv_cell_to_edges
1707 
1708 !----------------------------------------------------------------------
1709 
1710 subroutine r3_normalize(ax, ay, az)
1711 
1712  implicit none
1713 
1714  real(kind_real), intent(inout) :: ax, ay, az
1715  real(kind_real) :: mi
1716 
1717  mi = mpas_jedi_one_kr / sqrt(ax**2 + ay**2 + az**2)
1718 
1719  ax = ax * mi
1720  ay = ay * mi
1721  az = az * mi
1722 
1723 end subroutine r3_normalize
1724 
1725 !===============================================================================================================
1726 
1727 end module mpas4da_mod
1728 
subroutine, public da_self_mult(pool_a, zz)
Performs A = A * zz for pool A, zz a real number.
subroutine, public da_dot_product(pool_a, pool_b, dminfo, zprod)
Performs the dot_product given two pools of fields.
subroutine, public r3_normalize(ax, ay, az)
subroutine, public da_posdef(pool_a, fld_select)
Performs A = max(0.,A) for pool A.
subroutine, public da_constant(pool_a, realvalue, fld_select)
Performs A = constant. for pool A.
subroutine mpas_pool_template_field(srcFieldName, srcPool, dstFieldName, dstPool)
Add a field to dstPool that is templated on a srcPool field.
pure logical function match_scalar(scalarName, fieldName)
Test for a form of water vapor or hydrometeor of interest.
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.
character(len=1024) message
Definition: mpas4da_mod.F90:67
subroutine, public cvt_oopsmpas_date(inString2, outString2, iconv)
subroutine, public da_template_pool(geom, templatePool, nf, fieldnames)
Subset a pool from fields described in geom.
subroutine, public da_operator_addition(pool_a, pool_b)
Performs A = A + B for pools A and B.
integer function da_common_vars(pool_a, fieldname)
subroutine, public da_fldrms(pool_a, dminfo, fldrms, fld_select)
Performs basic statistics min/max/norm given a pool.
subroutine, public da_copy_all2sub_fields(domain, pool_a)
Performs a copy of allfield to a sub pool A.
subroutine, public uv_cell_to_edges(domain, u_field, v_field, du, lonCell, latCell, nCells, edgeNormalVectors, nEdgesOnCell, edgesOnCell, nVertLevels)
subroutine mpas_pool_demo(block)
Demonstrate basic usage of MPAS pools.
subroutine, public da_copy_sub2all_fields(domain, pool_a)
Performs a copy of a sub pool A to allfields.
subroutine, public da_setval(pool_a, zz)
Performs A = Val_R. for pool A.
pure logical function field_is_scalar(fieldName)
Test for a form of water vapor or hydrometeor of interest.
Definition: mpas4da_mod.F90:83
subroutine, public da_axpy(pool_a, pool_b, zz, fld_select)
Performs A = A + B * zz for pools A and B.
subroutine, public da_gpnorm(pool_a, dminfo, nf, pstat, fld_select)
Performs basic statistics min/max/norm given a pool.
subroutine, public da_random(pool_a, fld_select)
Performs random for pool A.
real(kind=kind_real), parameter mpas_jedi_half_kr
real(kind=kind_real), parameter mpas_jedi_zero_kr
real(kind=kind_real), parameter mpas_jedi_one_kr
integer function, dimension(:), allocatable, public getsolvedimsizes(pool, key, dimPool)
Returns an array with the dimension sizes of a pool field member.
logical function, public pool_has_field(pool, field)
Check for presence of field in pool.
Fortran derived type to hold geometry definition.