MPAS-JEDI
mpasjedi_linvarcha_c2a_mod.F90
Go to the documentation of this file.
1 ! (C) Copyright 2018 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 fckit_configuration_module, only: fckit_configuration
9 use fckit_log_module, only : fckit_log
10 use iso_c_binding
11 use kinds
12 
13 !oops
14 use oops_variables_mod, only: oops_variables
15 
16 !mpas-jedi
19 use mpas_geom_mod, only: mpas_geom
20 
21 !MPAS-Model
22 use mpas_derived_types
23 use mpas_field_routines
24 use mpas_pool_routines
25 use mpas_dmpar
26 
27 implicit none
28 
29 private
30 public :: mpasjedi_linvarcha_c2a
31 
32 !> Fortran derived type to hold configuration data for the B mat variable change
34  type (mpas_pool_type), pointer :: trajectories => null()
35  contains
36  procedure, public :: create
37  procedure, public :: delete
38  procedure, public :: multiply
39  procedure, public :: multiplyadjoint
40  procedure, public :: multiplyinverse
41  procedure, public :: multiplyinverseadjoint
43 
44 character(len=1024) :: message
45 
46 ! ------------------------------------------------------------------------------
47 contains
48 ! ------------------------------------------------------------------------------
49 
50 subroutine create(self, geom, bg, fg, f_conf, vars)
51 
52  implicit none
53  class(mpasjedi_linvarcha_c2a),intent(inout) :: self
54  type(mpas_fields), target, intent(in) :: bg
55  type(mpas_fields), target, intent(in) :: fg
56  type(mpas_geom), intent(in) :: geom
57  type(fckit_configuration), intent(in) :: f_conf
58  type(oops_variables), intent(in) :: vars
59 
60  type (field2DReal), pointer :: fld2d_t, fld2d_p, fld2d_qs
61 
62  integer :: ngrid
63 
64  if ( vars % has ('relhum') ) then
65  !-- set trajectories for linear variable change
66  ngrid = geom % nCells ! local + halo
67 
68  call mpas_pool_create_pool(self % trajectories)
69 
70  if ( .not. fg % has ('temperature') .or. .not. fg % has ('pressure') ) &
71  call abor1_ftn("LinVarChaC2A::LinVarChaC2A, mpasjedi_linvarcha_c2a_create: Trajectory failed")
72 
73  call fg%get('temperature', fld2d_t)
74  call fg%get('pressure' , fld2d_p)
75 
76  call mpas_duplicate_field(fld2d_t, fld2d_qs) ! for saturation specific humidity
77 
78  call da_tp_to_qs( fld2d_t % array(:,1:ngrid), fld2d_p % array(:,1:ngrid), fld2d_qs% array(:,1:ngrid))
79 
80  call mpas_pool_add_field(self % trajectories, 'spechum_sat', fld2d_qs)
81  end if
82 
83 end subroutine create
84 
85 ! ------------------------------------------------------------------------------
86 
87 subroutine delete(self)
88 
89  implicit none
90  class(mpasjedi_linvarcha_c2a),intent(inout) :: self
91 
92  if (associated(self % trajectories)) then
93  call mpas_pool_destroy_pool(self % trajectories)
94  end if
95 
96 end subroutine delete
97 
98 ! ------------------------------------------------------------------------------
99 
100 subroutine multiply(self,geom,xctl,xana)
101 
102  implicit none
103  class(mpasjedi_linvarcha_c2a),intent(inout) :: self
104  type(mpas_geom), intent(inout) :: geom
105  type(mpas_fields), intent(inout) :: xctl
106  type(mpas_fields), intent(inout) :: xana
107 
108  type (mpas_pool_iterator_type) :: poolItr
109  type (field2DReal), pointer :: fld2d_sf, fld2d_vp, fld2d_uRz, fld2d_uRm, fld2d_traj_qs
110  type (field2DReal), pointer :: fld2d_v_src, fld2d_sf_v, fld2d_e_src, fld2d_u
111 
112  real(kind=kind_real), dimension(:), pointer :: ptrr1_ctl, ptrr1_ana
113  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_ctl, ptrr2_ana
114  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_sh, ptrr2_rh
115 
116  integer :: k, ngrid, nCells, nVertices, nEdges
117 
118  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiply: xana % fldnames(:) =",xana % fldnames(:)
119  call fckit_log%debug(message)
120  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiply: xctl % fldnames(:) =",xctl % fldnames(:)
121  call fckit_log%debug(message)
122 
123  ngrid = geom % nCellsSolve ! local
124 
125  ! no variable change
126  call mpas_pool_begin_iteration(xctl%subFields)
127  do while ( mpas_pool_get_next_member(xctl%subFields, poolitr) )
128  ! NOTE: only 1d or 2d real fields are handled.
129  if (poolitr % memberType == mpas_pool_field .and. poolitr % dataType == mpas_pool_real) then
130  if( xana%has(poolitr % memberName) ) then ! check if the variable names are the same.
131  if (poolitr % nDims == 1) then
132  call xctl%get(poolitr % memberName, ptrr1_ctl)
133  call xana%get(poolitr % memberName, ptrr1_ana)
134  ptrr1_ana(1:ngrid)=ptrr1_ctl(1:ngrid)
135  else if (poolitr % nDims == 2) then
136  call xctl%get(poolitr % memberName, ptrr2_ctl)
137  call xana%get(poolitr % memberName, ptrr2_ana)
138  ptrr2_ana(:,1:ngrid)=ptrr2_ctl(:,1:ngrid)
139  else
140  write(message,*) '--> mpasjedi_linvarcha_c2a::multiply: poolItr % nDims == ',poolitr % nDims,' not handled'
141  call abor1_ftn(message)
142  end if
143  end if
144  end if
145  end do !- end of pool iteration
146 
147  if( xctl%has('relhum') .and. xana%has('spechum') ) then
148  call xctl%get( 'relhum', ptrr2_rh)
149  call xana%get( 'spechum', ptrr2_sh)
150  call mpas_pool_get_field(self%trajectories, 'spechum_sat', fld2d_traj_qs)
151  call pseudorh_to_spechum( ptrr2_rh(:,1:ngrid), ptrr2_sh(:,1:ngrid), fld2d_traj_qs%array(:,1:ngrid) )
152  end if
153 
154  if( xctl%has('stream_function') .and. xctl%has('velocity_potential') .and. &
155  xana%has('uReconstructZonal') .and. xana%has('uReconstructMeridional') ) then
156  call xctl%get( 'stream_function', fld2d_sf)
157  call xctl%get( 'velocity_potential', fld2d_vp)
158  call xana%get( 'uReconstructZonal', fld2d_urz)
159  call xana%get('uReconstructMeridional', fld2d_urm)
160 
161  call mpas_dmpar_exch_halo_field(fld2d_sf)
162  call mpas_dmpar_exch_halo_field(fld2d_vp)
163 
164  ncells = geom % nCells ! local + halo
165  nvertices = geom % nVertices ! local + halo
166  nedges = geom % nEdges ! local + halo
167 
168  ! duplicate two temporary working spaces
169  call mpas_pool_get_field( geom % domain % blocklist % allFields, 'vorticity', fld2d_v_src)
170  call mpas_duplicate_field(fld2d_v_src, fld2d_sf_v)
171  fld2d_sf_v % fieldName = 'stream_function at vertices'
172  call mpas_pool_get_field( geom % domain % blocklist % allFields, 'u', fld2d_e_src)
173  call mpas_duplicate_field(fld2d_e_src, fld2d_u)
174  fld2d_u % fieldName = 'Horizontal normal velocity at edges'
175 
176  call psichi_to_uv_edge_step1( geom, fld2d_sf%array(:,1:ncells), fld2d_sf_v%array(:,1:nvertices) )
177  call mpas_dmpar_exch_halo_field(fld2d_sf_v)
178  call psichi_to_uv_edge_step2( geom, fld2d_sf_v%array(:,1:nvertices), fld2d_vp%array(:,1:ncells), &
179  fld2d_u%array(:,1:nedges) )
180  call mpas_dmpar_exch_halo_field(fld2d_u)
181  call psichi_to_uv_edge_step3( geom, fld2d_u%array(:,1:nedges), &
182  fld2d_urz%array(:,1:ncells), fld2d_urm%array(:,1:ncells) )
183 
184  call mpas_deallocate_field(fld2d_sf_v)
185  call mpas_deallocate_field(fld2d_u)
186  end if
187 
188 end subroutine multiply
189 
190 ! ------------------------------------------------------------------------------
191 
192 subroutine multiplyadjoint(self,geom,xana,xctl)
193 
194  implicit none
195  class(mpasjedi_linvarcha_c2a),intent(inout) :: self
196  type(mpas_geom), intent(inout) :: geom
197  type(mpas_fields), intent(inout) :: xana
198  type(mpas_fields), intent(inout) :: xctl
199 
200  type (mpas_pool_iterator_type) :: poolItr
201  type (field2DReal), pointer :: fld2d_sf, fld2d_vp, fld2d_uRz, fld2d_uRm, fld2d_traj_qs
202  type (field2DReal), pointer :: fld2d_v_src, fld2d_sf_v, fld2d_e_src, fld2d_u
203 
204  real(kind=kind_real), dimension(:), pointer :: ptrr1_ctl, ptrr1_ana
205  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_ctl, ptrr2_ana
206  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_sh, ptrr2_rh
207  integer :: k, ngrid, nCells, nVertices, nEdges
208 
209  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiplyadjoint: xana % fldnames(:) =",xana % fldnames(:)
210  call fckit_log%debug(message)
211  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiplyadjoint: xctl % fldnames(:) =",xctl % fldnames(:)
212  call fckit_log%debug(message)
213 
214  ngrid = geom % nCellsSolve
215 
216  ! no variable change
217  call mpas_pool_begin_iteration(xctl%subFields)
218  do while ( mpas_pool_get_next_member(xctl%subFields, poolitr) )
219  ! NOTE: only 1d or 2d real fields are handled.
220  if (poolitr % memberType == mpas_pool_field .and. poolitr % dataType == mpas_pool_real) then
221  if( xana%has(poolitr % memberName) ) then ! check if the variable names are the same.
222  if (poolitr % nDims == 1) then
223  call xctl%get(poolitr % memberName, ptrr1_ctl)
224  call xana%get(poolitr % memberName, ptrr1_ana)
225  ptrr1_ctl(1:ngrid)=ptrr1_ctl(1:ngrid)+ptrr1_ana(1:ngrid)
226  ptrr1_ana(1:ngrid)=mpas_jedi_zero_kr
227  else if (poolitr % nDims == 2) then
228  call xctl%get(poolitr % memberName, ptrr2_ctl)
229  call xana%get(poolitr % memberName, ptrr2_ana)
230  ptrr2_ctl(:,1:ngrid)=ptrr2_ctl(:,1:ngrid)+ptrr2_ana(:,1:ngrid)
231  ptrr2_ana(:,1:ngrid)=mpas_jedi_zero_kr
232  else
233  write(message,*) '--> mpasjedi_linvarcha_c2a::multiplyadjoint: poolItr % nDims == ',poolitr % nDims,' not handled'
234  call abor1_ftn(message)
235  end if
236  end if
237  end if
238  end do !- end of pool iteration
239 
240  if( xctl%has('relhum') .and. xana%has('spechum') ) then
241  call xctl%get( 'relhum', ptrr2_rh)
242  call xana%get( 'spechum', ptrr2_sh)
243  call mpas_pool_get_field(self%trajectories, 'spechum_sat', fld2d_traj_qs)
244  call pseudorh_to_spechumad( ptrr2_rh(:,1:ngrid), ptrr2_sh(:,1:ngrid), fld2d_traj_qs%array(:,1:ngrid) )
245  end if
246 
247  if( xctl%has('stream_function') .and. xctl%has('velocity_potential') .and. &
248  xana%has('uReconstructZonal') .and. xana%has('uReconstructMeridional') ) then
249  call xctl%get( 'stream_function', fld2d_sf)
250  call xctl%get( 'velocity_potential', fld2d_vp)
251  call xana%get( 'uReconstructZonal', fld2d_urz)
252  call xana%get('uReconstructMeridional', fld2d_urm)
253 
254  ncells = geom % nCells ! local + halo
255  nvertices = geom % nVertices ! local + halo
256  nedges = geom % nEdges ! local + halo
257 
258  ! duplicate two temporary working spaces
259  call mpas_pool_get_field( geom % domain % blocklist % allFields, 'vorticity', fld2d_v_src)
260  call mpas_duplicate_field(fld2d_v_src, fld2d_sf_v)
261  fld2d_sf_v % fieldName = 'stream_function at vertices'
262  call mpas_pool_get_field( geom % domain % blocklist % allFields, 'u', fld2d_e_src)
263  call mpas_duplicate_field(fld2d_e_src, fld2d_u)
264  fld2d_u % fieldName = 'Horizontal normal velocity at edges'
265 
266  call psichi_to_uv_edge_step3ad( geom, fld2d_u%array(:,1:nedges), &
267  fld2d_urz%array(:,1:ncells), fld2d_urm%array(:,1:ncells) )
268  call mpas_dmpar_exch_halo_adj_field(fld2d_u)
269  call psichi_to_uv_edge_step2ad( geom, fld2d_sf_v%array(:,1:nvertices), fld2d_vp%array(:,1:ncells), &
270  fld2d_u%array(:,1:nedges) )
271  call mpas_dmpar_exch_halo_adj_field(fld2d_sf_v)
272  call psichi_to_uv_edge_step1ad( geom, fld2d_sf%array(:,1:ncells), fld2d_sf_v%array(:,1:nvertices) )
273 
274  call mpas_deallocate_field(fld2d_sf_v)
275  call mpas_deallocate_field(fld2d_u)
276 
277  fld2d_urz%array(:,1:geom%nCellsSolve)=mpas_jedi_zero_kr
278  fld2d_urm%array(:,1:geom%nCellsSolve)=mpas_jedi_zero_kr
279 
280  call mpas_dmpar_exch_halo_adj_field(fld2d_sf)
281  call mpas_dmpar_exch_halo_adj_field(fld2d_vp)
282  end if
283 
284 end subroutine multiplyadjoint
285 
286 ! ------------------------------------------------------------------------------
287 
288 subroutine multiplyinverse(self,geom,xana,xctl)
289 
290  implicit none
291  class(mpasjedi_linvarcha_c2a),intent(inout) :: self
292  type(mpas_geom), intent(inout) :: geom
293  type(mpas_fields), intent(inout) :: xana
294  type(mpas_fields), intent(inout) :: xctl
295 
296  type (mpas_pool_iterator_type) :: poolItr
297  type (field2DReal), pointer :: fld2d_traj_qs
298  real(kind=kind_real), dimension(:), pointer :: ptrr1_ctl, ptrr1_ana
299  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_ctl, ptrr2_ana
300  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_sh, ptrr2_rh
301  integer :: ngrid
302 
303  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiplyinverse: xana % fldnames(:) =",xana % fldnames(:)
304  call fckit_log%debug(message)
305  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiplyinverse: xctl % fldnames(:) =",xctl % fldnames(:)
306  call fckit_log%debug(message)
307 
308  ngrid = geom % nCellsSolve ! local
309 
310  ! no variable change
311  call mpas_pool_begin_iteration(xctl%subFields)
312  do while ( mpas_pool_get_next_member(xctl%subFields, poolitr) )
313  ! only 1d or 2d real fields are handled.
314  if (poolitr % memberType == mpas_pool_field .and. poolitr % dataType == mpas_pool_real) then
315  if( xana%has(poolitr % memberName) ) then ! check if the variable names are the same.
316  if (poolitr % nDims == 1) then
317  call xctl%get(poolitr % memberName, ptrr1_ctl)
318  call xana%get(poolitr % memberName, ptrr1_ana)
319  ptrr1_ctl(1:ngrid)=ptrr1_ana(1:ngrid)
320  else if (poolitr % nDims == 2) then
321  call xctl%get(poolitr % memberName, ptrr2_ctl)
322  call xana%get(poolitr % memberName, ptrr2_ana)
323  ptrr2_ctl(:,1:ngrid)=ptrr2_ana(:,1:ngrid)
324  else
325  write(message,*) '--> mpasjedi_linvarcha_c2a::multiplyinverse: poolItr % nDims == ',poolitr % nDims,' not handled'
326  call abor1_ftn(message)
327  end if
328  end if
329  end if
330  end do !- end of pool iteration
331 
332  if( xctl%has('relhum') .and. xana%has('spechum') ) then
333  call xana%get( 'spechum', ptrr2_sh)
334  call xctl%get( 'relhum', ptrr2_rh)
335  call mpas_pool_get_field(self%trajectories, 'spechum_sat', fld2d_traj_qs)
336  call pseudorh_to_spechum_inverse( ptrr2_rh(:,1:ngrid), ptrr2_sh(:,1:ngrid), fld2d_traj_qs%array(:,1:ngrid) )
337  end if
338 
339  !-- dummy inverse operator
340  if( xctl%has('stream_function') .and. xctl%has('velocity_potential') .and. &
341  xana%has('uReconstructZonal') .and. xana%has('uReconstructMeridional') ) then
342  call xana%get('uReconstructZonal', ptrr2_ana)
343  call xctl%get( 'stream_function', ptrr2_ctl)
344  ptrr2_ctl(:,1:ngrid)=ptrr2_ana(:,1:ngrid)
345 
346  call xana%get('uReconstructMeridional', ptrr2_ana)
347  call xctl%get( 'velocity_potential', ptrr2_ctl)
348  ptrr2_ctl(:,1:ngrid)=ptrr2_ana(:,1:ngrid)
349  end if
350 
351 end subroutine multiplyinverse
352 
353 ! ------------------------------------------------------------------------------
354 
355 subroutine multiplyinverseadjoint(self,geom,xctl,xana)
356 
357  implicit none
358  class(mpasjedi_linvarcha_c2a),intent(inout) :: self
359  type(mpas_geom), intent(inout) :: geom
360  type(mpas_fields), intent(inout) :: xctl
361  type(mpas_fields), intent(inout) :: xana
362 
363  type (mpas_pool_iterator_type) :: poolItr
364  type (field2DReal), pointer :: fld2d_traj_qs
365  real(kind=kind_real), dimension(:), pointer :: ptrr1_ctl, ptrr1_ana
366  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_ctl, ptrr2_ana
367  real(kind=kind_real), dimension(:,:), pointer :: ptrr2_sh, ptrr2_rh, ptrr2_traj_qs
368  integer :: ngrid
369 
370  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiplyinverseadjoint: xana % fldnames(:) =",xana % fldnames(:)
371  call fckit_log%debug(message)
372  write(message,*) "DEBUG: mpasjedi_linvarcha_c2a_multiplyinverseadjoint: xctl % fldnames(:) =",xctl % fldnames(:)
373  call fckit_log%debug(message)
374 
375  ngrid = geom % nCellsSolve ! local
376 
377  ! no variable change
378  call mpas_pool_begin_iteration(xctl%subFields)
379  do while ( mpas_pool_get_next_member(xctl%subFields, poolitr) )
380  ! only 1d or 2d real fields are handled.
381  if (poolitr % memberType == mpas_pool_field .and. poolitr % dataType == mpas_pool_real) then
382  if( xana%has(poolitr % memberName) ) then ! check if the variable names are the same.
383  if (poolitr % nDims == 1) then
384  call xctl%get(poolitr % memberName, ptrr1_ctl)
385  call xana%get(poolitr % memberName, ptrr1_ana)
386  ptrr1_ana(1:ngrid)=ptrr1_ana(1:ngrid)+ptrr1_ctl(1:ngrid)
387  ptrr1_ctl(1:ngrid)=mpas_jedi_zero_kr
388  else if (poolitr % nDims == 2) then
389  call xctl%get(poolitr % memberName, ptrr2_ctl)
390  call xana%get(poolitr % memberName, ptrr2_ana)
391  ptrr2_ana(:,1:ngrid)=ptrr2_ana(:,1:ngrid)+ptrr2_ctl(:,1:ngrid)
392  ptrr2_ctl(:,1:ngrid)=mpas_jedi_zero_kr
393  else
394  write(message,*) '--> mpasjedi_linvarcha_c2a::multiplyinverse: poolItr % nDims == ',poolitr % nDims,' not handled'
395  call abor1_ftn(message)
396  end if
397  end if
398  end if
399  end do !- end of pool iteration
400 
401  if( xctl%has('relhum') .and. xana%has('spechum') ) then
402  call xana%get( 'spechum', ptrr2_sh)
403  call xctl%get( 'relhum', ptrr2_rh)
404  call mpas_pool_get_field(self%trajectories, 'spechum_sat', fld2d_traj_qs)
405  call pseudorh_to_spechum_inversead( ptrr2_rh(:,1:ngrid), ptrr2_sh(:,1:ngrid), fld2d_traj_qs%array(:,1:ngrid) )
406  end if
407 
408  !-- dummy inverseAD operator
409  if( xctl%has('stream_function') .and. xctl%has('velocity_potential') .and. &
410  xana%has('uReconstructZonal') .and. xana%has('uReconstructMeridional') ) then
411  call xana%get('uReconstructZonal', ptrr2_ana)
412  call xctl%get( 'stream_function', ptrr2_ctl)
413  ptrr2_ana(:,1:ngrid)=ptrr2_ana(:,1:ngrid)+ptrr2_ctl(:,1:ngrid)
414  ptrr2_ctl(:,1:ngrid)=mpas_jedi_zero_kr
415 
416  call xana%get('uReconstructMeridional', ptrr2_ana)
417  call xctl%get( 'velocity_potential', ptrr2_ctl)
418  ptrr2_ana(:,1:ngrid)=ptrr2_ana(:,1:ngrid)+ptrr2_ctl(:,1:ngrid)
419  ptrr2_ctl(:,1:ngrid)=mpas_jedi_zero_kr
420  end if
421 
422 end subroutine multiplyinverseadjoint
423 
424 !-------------------------------------------------------------------------------
425 ! Adjoint code of subroutine mpas_reconstruct()
426 ! in MPAS-Model/src/operators/mpas_vector_reconstruction.F
427 !
428 subroutine mpas_reconstruct_2dad(meshPool, u, uReconstructX, uReconstructY, uReconstructZ,&
429  uReconstructZonal, uReconstructMeridional, includeHalos)!{{{
430 
431  implicit none
432 
433  type (mpas_pool_type), intent(in) :: meshPool !< Input: Mesh information
434  real (kind=rkind), dimension(:,:), intent(inout) :: u !< Input: Velocity field on edges
435  real (kind=rkind), dimension(:,:), intent(inout) :: ureconstructx !< Output: X Component of velocity reconstructed to cell centers
436  real (kind=rkind), dimension(:,:), intent(inout) :: ureconstructy !< Output: Y Component of velocity reconstructed to cell centers
437  real (kind=rkind), dimension(:,:), intent(inout) :: ureconstructz !< Output: Z Component of velocity reconstructed to cell centers
438  real (kind=rkind), dimension(:,:), intent(inout) :: ureconstructzonal !< Output: Zonal Component of velocity reconstructed to cell centers
439  real (kind=rkind), dimension(:,:), intent(inout) :: ureconstructmeridional !< Output: Meridional Component of velocity reconstructed to cell centers
440  logical, optional, intent(in) :: includeHalos !< Input: Optional logical that allows reconstruction over halo regions
441 
442  ! temporary arrays needed in the compute procedure
443  logical :: includeHalosLocal
444  integer, pointer :: nCells
445  integer, dimension(:,:), pointer :: edgesOnCell
446  integer, dimension(:), pointer :: nEdgesOnCell
447  integer :: iCell,iEdge, i
448  real(kind=rkind), dimension(:), pointer :: latcell, loncell
449 
450  real (kind=rkind), dimension(:,:,:), pointer :: coeffs_reconstruct
451 
452  logical, pointer :: on_a_sphere
453 
454  real (kind=rkind) :: clat, slat, clon, slon
455 
456  if ( present(includehalos) ) then
457  includehaloslocal = includehalos
458  else
459  includehaloslocal = .false.
460  end if
461 
462  ! stored arrays used during compute procedure
463  call mpas_pool_get_array(meshpool, 'coeffs_reconstruct', coeffs_reconstruct)
464 
465  ! temporary variables
466  call mpas_pool_get_array(meshpool, 'nEdgesOnCell', nedgesoncell)
467  call mpas_pool_get_array(meshpool, 'edgesOnCell', edgesoncell)
468 
469  if ( includehaloslocal ) then
470  call mpas_pool_get_dimension(meshpool, 'nCells', ncells)
471  else
472  call mpas_pool_get_dimension(meshpool, 'nCellsSolve', ncells)
473  end if
474 
475  call mpas_pool_get_array(meshpool, 'latCell', latcell)
476  call mpas_pool_get_array(meshpool, 'lonCell', loncell)
477 
478  call mpas_pool_get_config(meshpool, 'on_a_sphere', on_a_sphere)
479 
480  ! initialize the reconstructed vectors
481  ureconstructx = mpas_jedi_zero_kr
482  ureconstructy = mpas_jedi_zero_kr
483  ureconstructz = mpas_jedi_zero_kr
484  if (on_a_sphere) then
485  do icell = 1, ncells
486  clat = cos(latcell(icell))
487  slat = sin(latcell(icell))
488  clon = cos(loncell(icell))
489  slon = sin(loncell(icell))
490 
491  ureconstructx(:,icell) = ureconstructx(:,icell) - clon*slat * ureconstructmeridional(:,icell)
492  ureconstructy(:,icell) = ureconstructy(:,icell) - slon*slat * ureconstructmeridional(:,icell)
493  ureconstructz(:,icell) = ureconstructz(:,icell) + clat * ureconstructmeridional(:,icell)
494 
495  ureconstructx(:,icell) = ureconstructx(:,icell) - slon * ureconstructzonal(:,icell)
496  ureconstructy(:,icell) = ureconstructy(:,icell) + clon * ureconstructzonal(:,icell)
497  end do
498  else
499  do icell = 1, ncells
500  ureconstructy(:,icell) = ureconstructy(:,icell) + ureconstructmeridional(:,icell)
501  ureconstructx(:,icell) = ureconstructx(:,icell) + ureconstructzonal(:,icell)
502  end do
503  end if
504 
505  ! loop over cell centers
506  do icell = 1, ncells
507  ! a more efficient reconstruction where rbf_values*matrix_reconstruct
508  ! has been precomputed in coeffs_reconstruct
509  do i=1,nedgesoncell(icell)
510  iedge = edgesoncell(i,icell)
511  u(:,iedge) = u(:,iedge) + coeffs_reconstruct(1,i,icell) * ureconstructx(:,icell) &
512  + coeffs_reconstruct(2,i,icell) * ureconstructy(:,icell) &
513  + coeffs_reconstruct(3,i,icell) * ureconstructz(:,icell)
514  enddo
515  enddo ! iCell
516 
517 end subroutine mpas_reconstruct_2dad!}}}
518 
519 !-------------------------------------------------------------------------------
520 
521 subroutine mpas_reconstruct_1dad(meshPool, u, uReconstructX, uReconstructY, uReconstructZ,&
522  uReconstructZonal, uReconstructMeridional, includeHalos)!{{{
523 
524  implicit none
525 
526  type (mpas_pool_type), intent(in) :: meshPool !< Input: Mesh information
527  real (kind=rkind), dimension(:), intent(inout) :: u !< Input: Velocity field on edges
528  real (kind=rkind), dimension(:), intent(inout) :: ureconstructx !< Output: X Component of velocity reconstructed to cell centers
529  real (kind=rkind), dimension(:), intent(inout) :: ureconstructy !< Output: Y Component of velocity reconstructed to cell centers
530  real (kind=rkind), dimension(:), intent(inout) :: ureconstructz !< Output: Z Component of velocity reconstructed to cell centers
531  real (kind=rkind), dimension(:), intent(inout) :: ureconstructzonal !< Output: Zonal Component of velocity reconstructed to cell centers
532  real (kind=rkind), dimension(:), intent(inout) :: ureconstructmeridional !< Output: Meridional Component of velocity reconstructed to cell centers
533  logical, optional, intent(in) :: includeHalos !< Input: Optional logical that allows reconstruction over halo regions
534 
535  ! temporary arrays needed in the compute procedure
536  logical :: includeHalosLocal
537  integer, pointer :: nCells
538  integer, dimension(:,:), pointer :: edgesOnCell
539  integer, dimension(:), pointer :: nEdgesOnCell
540  integer :: iCell,iEdge, i
541  real(kind=rkind), dimension(:), pointer :: latcell, loncell
542 
543  real (kind=rkind), dimension(:,:,:), pointer :: coeffs_reconstruct
544 
545  logical, pointer :: on_a_sphere
546 
547  real (kind=rkind) :: clat, slat, clon, slon
548 
549  if ( present(includehalos) ) then
550  includehaloslocal = includehalos
551  else
552  includehaloslocal = .false.
553  end if
554 
555  ! stored arrays used during compute procedure
556  call mpas_pool_get_array(meshpool, 'coeffs_reconstruct', coeffs_reconstruct)
557 
558  ! temporary variables
559  call mpas_pool_get_array(meshpool, 'nEdgesOnCell', nedgesoncell)
560  call mpas_pool_get_array(meshpool, 'edgesOnCell', edgesoncell)
561 
562  if ( includehaloslocal ) then
563  call mpas_pool_get_dimension(meshpool, 'nCells', ncells)
564  else
565  call mpas_pool_get_dimension(meshpool, 'nCellsSolve', ncells)
566  end if
567 
568  call mpas_pool_get_array(meshpool, 'latCell', latcell)
569  call mpas_pool_get_array(meshpool, 'lonCell', loncell)
570 
571  call mpas_pool_get_config(meshpool, 'on_a_sphere', on_a_sphere)
572 
573  ! initialize the reconstructed vectors
574  ureconstructx = mpas_jedi_zero_kr
575  ureconstructy = mpas_jedi_zero_kr
576  ureconstructz = mpas_jedi_zero_kr
577  if (on_a_sphere) then
578  do icell = 1, ncells
579  clat = cos(latcell(icell))
580  slat = sin(latcell(icell))
581  clon = cos(loncell(icell))
582  slon = sin(loncell(icell))
583 
584  ureconstructx(icell) = ureconstructx(icell) - clon*slat * ureconstructmeridional(icell)
585  ureconstructy(icell) = ureconstructy(icell) - slon*slat * ureconstructmeridional(icell)
586  ureconstructz(icell) = ureconstructz(icell) + clat * ureconstructmeridional(icell)
587 
588  ureconstructx(icell) = ureconstructx(icell) - slon * ureconstructzonal(icell)
589  ureconstructy(icell) = ureconstructy(icell) + clon * ureconstructzonal(icell)
590  end do
591  else
592  do icell = 1, ncells
593  ureconstructy(icell) = ureconstructy(icell) + ureconstructmeridional(icell)
594  ureconstructx(icell) = ureconstructx(icell) + ureconstructzonal(icell)
595  end do
596  end if
597 
598  ! loop over cell centers
599  do icell = 1, ncells
600  ! a more efficient reconstruction where rbf_values*matrix_reconstruct
601  ! has been precomputed in coeffs_reconstruct
602  do i=1,nedgesoncell(icell)
603  iedge = edgesoncell(i,icell)
604  u(iedge) = u(iedge) + coeffs_reconstruct(1,i,icell) * ureconstructx(icell) &
605  + coeffs_reconstruct(2,i,icell) * ureconstructy(icell) &
606  + coeffs_reconstruct(3,i,icell) * ureconstructz(icell)
607  enddo
608  enddo ! iCell
609 
610 end subroutine mpas_reconstruct_1dad!}}}
611 
612 !-------------------------------------------------------------------------------
613 ! input : psi & chi @ cell center
614 ! output : u & v @ cell center
615 subroutine psichi_to_uv_center(geom, psi, chi, u, v)
616 
617  implicit none
618  type (mpas_geom), intent(in) :: geom !< geometry
619  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(in) :: psi, chi
620  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(out) :: u, v
621 
622  integer :: iC, iE, j
623  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCellsSolve) :: &
624  psi_line_intg_dx, psi_line_intg_dy, chi_line_intg_dx, chi_line_intg_dy
625 
626  psi_line_intg_dx=mpas_jedi_zero_kr
627  psi_line_intg_dy=mpas_jedi_zero_kr
628  chi_line_intg_dx=mpas_jedi_zero_kr
629  chi_line_intg_dy=mpas_jedi_zero_kr
630 
631  do ic = 1, geom%nCellsSolve
632  do j = 1, geom%nEdgesOnCell(ic) ! or geom%maxEdges
633  ie = geom%edgesOnCell(j,ic)
634  psi_line_intg_dx(:,ic) = psi_line_intg_dx(:,ic) &
635  + mpas_jedi_half_kr * ( psi(:,geom%cellsOnEdge(1,ie)) + psi(:,geom%cellsOnEdge(2,ie)) ) &
636  * geom%dvEdge(ie) * sin(mpas_jedi_zero_kr-geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic)
637  psi_line_intg_dy(:,ic) = psi_line_intg_dy(:,ic) &
638  + mpas_jedi_half_kr * ( psi(:,geom%cellsOnEdge(1,ie)) + psi(:,geom%cellsOnEdge(2,ie)) ) &
639  * geom%dvEdge(ie) * cos(geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic)
640  chi_line_intg_dx(:,ic) = chi_line_intg_dx(:,ic) &
641  + mpas_jedi_half_kr * ( chi(:,geom%cellsOnEdge(1,ie)) + chi(:,geom%cellsOnEdge(2,ie)) ) &
642  * geom%dvEdge(ie) * sin(mpas_jedi_zero_kr-geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic)
643  chi_line_intg_dy(:,ic) = chi_line_intg_dy(:,ic) &
644  + mpas_jedi_half_kr * ( chi(:,geom%cellsOnEdge(1,ie)) + chi(:,geom%cellsOnEdge(2,ie)) ) &
645  * geom%dvEdge(ie) * cos(geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic)
646  enddo !- j
647  enddo !- iC
648 
649  do ic=1, geom%nCellsSolve
650  u(:,ic) = ( psi_line_intg_dx(:,ic) - chi_line_intg_dy(:,ic) ) / geom%areaCell(ic)
651  v(:,ic) = ( psi_line_intg_dy(:,ic) + chi_line_intg_dx(:,ic) ) / geom%areaCell(ic)
652  enddo
653 
654 end subroutine psichi_to_uv_center
655 
656 !-------------------------------------------------------------------------------
657 
658 subroutine psichi_to_uv_centerad(geom, psi, chi, u, v)
659 
660  implicit none
661  type (mpas_geom), intent(in) :: geom !< geometry
662  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(inout) :: psi, chi
663  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(inout) :: u, v
664 
665  integer :: iC, iE, j
666  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCellsSolve) :: &
667  psi_line_intg_dx, psi_line_intg_dy, chi_line_intg_dx, chi_line_intg_dy
668 
669  psi_line_intg_dx=mpas_jedi_zero_kr
670  psi_line_intg_dy=mpas_jedi_zero_kr
671  chi_line_intg_dx=mpas_jedi_zero_kr
672  chi_line_intg_dy=mpas_jedi_zero_kr
673 
674  do ic=1, geom%nCellsSolve
675  psi_line_intg_dx(:,ic) = psi_line_intg_dx(:,ic) + u(:,ic) / geom%areaCell(ic)
676  chi_line_intg_dy(:,ic) = chi_line_intg_dy(:,ic) - u(:,ic) / geom%areaCell(ic)
677  psi_line_intg_dy(:,ic) = psi_line_intg_dy(:,ic) + v(:,ic) / geom%areaCell(ic)
678  chi_line_intg_dx(:,ic) = chi_line_intg_dx(:,ic) + v(:,ic) / geom%areaCell(ic)
679  enddo
680 
681  do ic = 1, geom%nCellsSolve
682  do j = 1, geom%nEdgesOnCell(ic) ! or geom%maxEdges
683  ie = geom%edgesOnCell(j,ic)
684 
685  psi(:,geom%cellsOnEdge(1,ie)) = psi(:,geom%cellsOnEdge(1,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
686  * sin(mpas_jedi_zero_kr-geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * psi_line_intg_dx(:,ic)
687  psi(:,geom%cellsOnEdge(2,ie)) = psi(:,geom%cellsOnEdge(2,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
688  * sin(mpas_jedi_zero_kr-geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * psi_line_intg_dx(:,ic)
689 
690  psi(:,geom%cellsOnEdge(1,ie)) = psi(:,geom%cellsOnEdge(1,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
691  * cos(geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * psi_line_intg_dy(:,ic)
692  psi(:,geom%cellsOnEdge(2,ie)) = psi(:,geom%cellsOnEdge(2,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
693  * cos(geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * psi_line_intg_dy(:,ic)
694 
695  chi(:,geom%cellsOnEdge(1,ie)) = chi(:,geom%cellsOnEdge(1,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
696  * sin(mpas_jedi_zero_kr-geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * chi_line_intg_dx(:,ic)
697  chi(:,geom%cellsOnEdge(2,ie)) = chi(:,geom%cellsOnEdge(2,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
698  * sin(mpas_jedi_zero_kr-geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * chi_line_intg_dx(:,ic)
699 
700  chi(:,geom%cellsOnEdge(1,ie)) = chi(:,geom%cellsOnEdge(1,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
701  * cos(geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * chi_line_intg_dy(:,ic)
702  chi(:,geom%cellsOnEdge(2,ie)) = chi(:,geom%cellsOnEdge(2,ie)) + mpas_jedi_half_kr * geom%dvEdge(ie)&
703  * cos(geom%angleEdge(ie)) * geom%edgesOnCell_sign(j,ic) * chi_line_intg_dy(:,ic)
704  enddo !- j
705  enddo !- iC
706 
707 end subroutine psichi_to_uv_centerad
708 
709 !-------------------------------------------------------------------------------
710 ! input : psi @ cell center
711 ! output : psi @ vertices
712 subroutine psichi_to_uv_edge_step1(geom, psi, psi_v)
713 
714  implicit none
715  type (mpas_geom), intent(in) :: geom !< geometry
716  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(in) :: psi
717  real (kind=kind_real), dimension(geom%nVertLevels,geom%nVertices), intent(inout) :: psi_v
718 
719  integer :: iC, iV, j
720 
721  psi_v=mpas_jedi_zero_kr
722 
723  ! Interpolate psi in cell center to vertice
724  do iv = 1, geom%nVerticesSolve ! local
725  do j = 1, geom%vertexDegree
726  ic = geom%cellsOnVertex(j,iv)
727  psi_v(:,iv) = psi_v(:,iv) + geom%kiteAreasOnVertex(j,iv) * psi(:,ic)
728  enddo
729  psi_v(:,iv) = psi_v(:,iv) / geom%areaTriangle(iv)
730  enddo
731 
732 end subroutine psichi_to_uv_edge_step1
733 
734 !-------------------------------------------------------------------------------
735 ! input : psi @ vertices, chi @ cell center
736 ! output : edge_normal_wind @ edges
737 subroutine psichi_to_uv_edge_step2(geom, psi_v, chi, edge_normal_wind)
738 
739  implicit none
740  type (mpas_geom), intent(in) :: geom !< geometry
741  real (kind=kind_real), dimension(geom%nVertLevels,geom%nVertices), intent(in) :: psi_v
742  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(in) :: chi
743  real (kind=kind_real), dimension(geom%nVertLevels,geom%nEdges), intent(inout) :: edge_normal_wind
744 
745  integer :: iE
746 
747  edge_normal_wind=mpas_jedi_zero_kr
748 
749  !get edge_normal_wind
750  do ie = 1, geom%nEdgesSolve ! local
751  edge_normal_wind(:,ie) = edge_normal_wind(:,ie) - &
752  ( chi(:,geom%cellsOnEdge(2,ie)) - chi(:,geom%cellsOnEdge(1,ie)) ) / geom%dcEdge(ie) - &
753  ( psi_v(:,geom%verticesOnEdge(2,ie)) - psi_v(:,geom%verticesOnEdge(1,ie)) ) / geom%dvEdge(ie)
754  enddo
755 
756 end subroutine psichi_to_uv_edge_step2
757 
758 !-------------------------------------------------------------------------------
759 ! input : edge_normal_wind @ edges
760 ! output : zonal- and meridional- wind @ cell center
761 subroutine psichi_to_uv_edge_step3(geom, edge_normal_wind, u, v)
762 
763  use mpas_vector_reconstruction
764 
765  implicit none
766  type (mpas_geom), intent(in) :: geom !< geometry
767  real (kind=kind_real), dimension(geom%nVertLevels,geom%nEdges), intent(in) :: edge_normal_wind
768  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(inout) :: u, v
769 
770  real (kind=kind_real), dimension(:,:), allocatable :: &
771  ureconstructx, ureconstructy, ureconstructz
772  type (mpas_pool_type), pointer :: mesh
773 
774  allocate(ureconstructx(geom%nVertLevels,geom%nCells))
775  allocate(ureconstructy(geom%nVertLevels,geom%nCells))
776  allocate(ureconstructz(geom%nVertLevels,geom%nCells))
777 
778  !edge_wind -> u/v @ cell center
779  call mpas_pool_get_subpool( geom % domain % blocklist % structs, 'mesh', mesh)
780  call mpas_reconstruct(mesh, edge_normal_wind, &
781  ureconstructx, &
782  ureconstructy, &
783  ureconstructz, &
784  u, &
785  v, .false. ) ! local only, no halo calculation
786 
787  deallocate(ureconstructx, ureconstructy, ureconstructz)
788 
789 end subroutine psichi_to_uv_edge_step3
790 
791 !-------------------------------------------------------------------------------
792 subroutine psichi_to_uv_edge_step1ad(geom, psi, psi_v)
793 
794  implicit none
795  type (mpas_geom), intent(in) :: geom !< geometry
796  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(inout) :: psi
797  real (kind=kind_real), dimension(geom%nVertLevels,geom%nVertices), intent(inout) :: psi_v
798 
799  integer :: iC, iV, j
800 
802 
803  ! Interpolate psi in cell center to vertice
804  do iv = 1, geom%nVerticesSolve ! local
805  psi_v(:,iv) = psi_v(:,iv) / geom%areaTriangle(iv)
806  do j = 1, geom%vertexDegree
807  ic = geom%cellsOnVertex(j,iv)
808  psi(:,ic) = psi(:,ic) + geom%kiteAreasOnVertex(j,iv) * psi_v(:,iv)
809  enddo
810  enddo
811 
812 end subroutine psichi_to_uv_edge_step1ad
813 
814 !-------------------------------------------------------------------------------
815 subroutine psichi_to_uv_edge_step2ad(geom, psi_v, chi, edge_normal_wind)
816 
817  implicit none
818  type (mpas_geom), intent(in) :: geom !< geometry
819  real (kind=kind_real), dimension(geom%nVertLevels,geom%nVertices), intent(inout) :: psi_v
820  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(inout) :: chi
821  real (kind=kind_real), dimension(geom%nVertLevels,geom%nEdges), intent(inout) :: edge_normal_wind
822 
823  integer :: iE
824 
826  psi_v=mpas_jedi_zero_kr
827 
828  !get edge_normal_wind
829  do ie = 1, geom%nEdgesSolve ! local
830  chi(:,geom%cellsOnEdge(2,ie)) = chi(:,geom%cellsOnEdge(2,ie)) - edge_normal_wind(:,ie) / geom%dcEdge(ie)
831  chi(:,geom%cellsOnEdge(1,ie)) = chi(:,geom%cellsOnEdge(1,ie)) + edge_normal_wind(:,ie) / geom%dcEdge(ie)
832  psi_v(:,geom%verticesOnEdge(2,ie)) = psi_v(:,geom%verticesOnEdge(2,ie)) - edge_normal_wind(:,ie) / geom%dvEdge(ie)
833  psi_v(:,geom%verticesOnEdge(1,ie)) = psi_v(:,geom%verticesOnEdge(1,ie)) + edge_normal_wind(:,ie) / geom%dvEdge(ie)
834  enddo
835 
836 end subroutine psichi_to_uv_edge_step2ad
837 
838 !-------------------------------------------------------------------------------
839 subroutine psichi_to_uv_edge_step3ad(geom, edge_normal_wind, u, v)
840 
841  use mpas_vector_reconstruction
842 
843  implicit none
844  type (mpas_geom), intent(in) :: geom !< geometry
845  real (kind=kind_real), dimension(geom%nVertLevels,geom%nEdges), intent(inout) :: edge_normal_wind
846  real (kind=kind_real), dimension(geom%nVertLevels,geom%nCells), intent(inout) :: u, v
847 
848  real (kind=kind_real), dimension(:,:), allocatable :: &
849  ureconstructx, ureconstructy, ureconstructz
850  type (mpas_pool_type), pointer :: mesh
851 
852  edge_normal_wind=mpas_jedi_zero_kr
853 
854  allocate(ureconstructx(geom%nVertLevels,geom%nCells))
855  allocate(ureconstructy(geom%nVertLevels,geom%nCells))
856  allocate(ureconstructz(geom%nVertLevels,geom%nCells))
857 
858  call mpas_pool_get_subpool( geom % domain % blocklist % structs, 'mesh', mesh)
859  call mpas_reconstruct_2dad(mesh, edge_normal_wind, &
860  ureconstructx, &
861  ureconstructy, &
862  ureconstructz, &
863  u, &
864  v, .false. ) ! local only, no halo calculation
865 
866  deallocate(ureconstructx, ureconstructy, ureconstructz)
867 
868 end subroutine psichi_to_uv_edge_step3ad
869 
870 ! ------------------------------------------------------------------------------
871 
872 elemental subroutine pseudorh_to_spechum(pseudorh,spechum,saturation_spechum)
873 
874  implicit none
875  real (kind=kind_real), intent(in) :: pseudorh
876  real (kind=kind_real), intent(out) :: spechum
877  real (kind=kind_real), intent(in) :: saturation_spechum
878 
879  spechum = pseudorh * saturation_spechum / mpas_jedi_hundred_kr
880 
881 end subroutine pseudorh_to_spechum
882 
883 elemental subroutine pseudorh_to_spechumad(pseudorh,spechum,saturation_spechum)
884 
885  implicit none
886  real (kind=kind_real), intent(inout) :: pseudorh
887  real (kind=kind_real), intent(inout) :: spechum
888  real (kind=kind_real), intent(in) :: saturation_spechum
889 
890  pseudorh = pseudorh + spechum * saturation_spechum / mpas_jedi_hundred_kr
891  spechum = mpas_jedi_zero_kr
892 
893 end subroutine pseudorh_to_spechumad
894 
895 elemental subroutine pseudorh_to_spechum_inverse(pseudorh,spechum,saturation_spechum)
896 
897  implicit none
898  real (kind=kind_real), intent(out) :: pseudorh
899  real (kind=kind_real), intent(in) :: spechum
900  real (kind=kind_real), intent(in) :: saturation_spechum
901 
902  pseudorh = spechum / saturation_spechum * mpas_jedi_hundred_kr
903 
904 end subroutine pseudorh_to_spechum_inverse
905 
906 elemental subroutine pseudorh_to_spechum_inversead(pseudorh,spechum,saturation_spechum)
907 
908  implicit none
909  real (kind=kind_real), intent(inout) :: pseudorh
910  real (kind=kind_real), intent(inout) :: spechum
911  real (kind=kind_real), intent(in) :: saturation_spechum
912 
913  spechum = spechum + pseudorh / saturation_spechum * mpas_jedi_hundred_kr
914  pseudorh = mpas_jedi_zero_kr
915 
916 end subroutine pseudorh_to_spechum_inversead
917 
918 elemental subroutine da_tp_to_qs( t, p, qs)
919 
920  !---------------------------------------------------------------------------
921  ! Purpose: Convert T/p to saturation specific humidity.
922  !
923  ! Method: qs = es_alpha * es / ( p - ( 1 - rd_over_rv ) * es ).
924  ! use Rogers & Yau (1989) formula: es = a exp( bTc / (T_c + c) )
925  !--------------------------------------------------------------------------
926 
927  implicit none
928 
929  real (kind=kind_real), intent(in) :: t, p
930  real (kind=kind_real), intent(out) :: qs
931 
932  real (kind=kind_real) :: es
933  real (kind=kind_real) :: t_c ! T in degreesC.
934 
935  !---------------------------------------------------------------------------
936  ! [1.0] initialise:
937  !---------------------------------------------------------------------------
938  t_c = t - t_kelvin
939 
940  !---------------------------------------------------------------------------
941  ! [2.0] Calculate saturation vapour pressure:
942  !---------------------------------------------------------------------------
943  es = es_alpha * exp( es_beta * t_c / ( t_c + es_gamma ) )
944 
945  !---------------------------------------------------------------------------
946  ! [3.0] Calculate saturation specific humidity:
947  !---------------------------------------------------------------------------
948  qs = rd_over_rv * es / ( p - rd_over_rv1 * es )
949 
950 end subroutine da_tp_to_qs
951 
real(kind=kind_real), parameter mpas_jedi_half_kr
real(kind=kind_real), parameter rd_over_rv1
real(kind=kind_real), parameter t_kelvin
real(kind=kind_real), parameter es_gamma
real(kind=kind_real), parameter es_alpha
real(kind=kind_real), parameter mpas_jedi_zero_kr
real(kind=kind_real), parameter es_beta
real(kind=kind_real), parameter rd_over_rv
real(kind=kind_real), parameter mpas_jedi_hundred_kr
subroutine mpas_reconstruct_1dad(meshPool, u, uReconstructX, uReconstructY, uReconstructZ, uReconstructZonal, uReconstructMeridional, includeHalos)
elemental subroutine pseudorh_to_spechumad(pseudorh, spechum, saturation_spechum)
subroutine psichi_to_uv_edge_step2ad(geom, psi_v, chi, edge_normal_wind)
subroutine psichi_to_uv_edge_step1(geom, psi, psi_v)
subroutine create(self, geom, bg, fg, f_conf, vars)
subroutine multiply(self, geom, xctl, xana)
subroutine multiplyadjoint(self, geom, xana, xctl)
subroutine psichi_to_uv_edge_step2(geom, psi_v, chi, edge_normal_wind)
subroutine multiplyinverse(self, geom, xana, xctl)
subroutine mpas_reconstruct_2dad(meshPool, u, uReconstructX, uReconstructY, uReconstructZ, uReconstructZonal, uReconstructMeridional, includeHalos)
elemental subroutine pseudorh_to_spechum_inverse(pseudorh, spechum, saturation_spechum)
elemental subroutine da_tp_to_qs(t, p, qs)
subroutine psichi_to_uv_edge_step1ad(geom, psi, psi_v)
subroutine psichi_to_uv_centerad(geom, psi, chi, u, v)
subroutine psichi_to_uv_center(geom, psi, chi, u, v)
subroutine psichi_to_uv_edge_step3ad(geom, edge_normal_wind, u, v)
elemental subroutine pseudorh_to_spechum(pseudorh, spechum, saturation_spechum)
elemental subroutine pseudorh_to_spechum_inversead(pseudorh, spechum, saturation_spechum)
subroutine multiplyinverseadjoint(self, geom, xctl, xana)
subroutine psichi_to_uv_edge_step3(geom, edge_normal_wind, u, v)
Fortran derived type to hold MPAS field.
Fortran derived type to hold geometry definition.
Fortran derived type to hold configuration data for the B mat variable change.