FV3-JEDI
fv3jedi_linvarcha_a2m_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2018-2019 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 iso_c_binding
9 use fckit_configuration_module, only: fckit_configuration
11 
15 
17 
18 use wind_vt_mod, only: a2d, d2a, d2a_ad, a2d_ad
19 
20 implicit none
21 private
22 
23 public :: fv3jedi_linvarcha_a2m
24 public :: create
25 public :: delete
26 public :: multiply
27 public :: multiplyadjoint
28 public :: multiplyinverse
29 public :: multiplyinverseadjoint
30 
32  integer :: dummy
33 end type fv3jedi_linvarcha_a2m
34 
35 ! ------------------------------------------------------------------------------
36 
37 contains
38 
39 ! ------------------------------------------------------------------------------
40 
41 subroutine create(self, geom, bg, fg, c_conf)
42 
43 implicit none
44 type(fv3jedi_linvarcha_a2m), intent(inout) :: self
45 type(fv3jedi_geom), intent(inout) :: geom
46 type(fv3jedi_state), intent(in) :: bg
47 type(fv3jedi_state), intent(in) :: fg
48 type(c_ptr), intent(in) :: c_conf
49 
50 type(fckit_configuration) :: f_conf
51 
52 f_conf = fckit_configuration(c_conf)
53 
54 !Nothing yet as transforms are linear
55 
56 end subroutine create
57 
58 ! ------------------------------------------------------------------------------
59 
60 subroutine delete(self)
61 
62 implicit none
63 type(fv3jedi_linvarcha_a2m), intent(inout) :: self
64 
65 end subroutine delete
66 
67 ! ------------------------------------------------------------------------------
68 
69 subroutine multiply(self,geom,xana,xmod)
70 
71 implicit none
72 type(fv3jedi_linvarcha_a2m), intent(inout) :: self
73 type(fv3jedi_geom), intent(inout) :: geom
74 type(fv3jedi_increment), intent(in) :: xana
75 type(fv3jedi_increment), intent(inout) :: xmod
76 
77 integer :: index_ana, index_mod, index_ana_found
78 integer :: k
79 logical :: failed
80 
81 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ua
82 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_va
83 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ps
84 
85 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_ud
86 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_vd
87 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_delp
88 
89 do index_mod = 1, xmod%nf
90 
91  index_ana_found = -1
92  failed = .true.
93 
94  !Check analysis for presence of field
95  do index_ana = 1, xana%nf
96  if (xmod%fields(index_mod)%fv3jedi_name == xana%fields(index_ana)%fv3jedi_name) then
97  index_ana_found = index_ana
98  exit
99  endif
100  enddo
101 
102  if (index_ana_found >= 0) then
103 
104  !OK, direct copy
105  xmod%fields(index_mod)%array = xana%fields(index_ana_found)%array
106  failed = .false.
107  if (xmod%f_comm%rank() == 0) write(*,"(A, A10, A, A10)") &
108  "A2M Multiply: analysis increment "//xana%fields(index_ana_found)%fv3jedi_name(1:10)&
109  //" => linearized model "//xmod%fields(index_mod)%fv3jedi_name(1:10)
110 
111  elseif (xmod%fields(index_mod)%fv3jedi_name == 'ud') then
112 
113  !Special case: A-grid analysis, D-Grid model
114  if (xana%has_field('ua')) then
115  call xana%get_field('ua', xana_ua)
116  call xana%get_field('va', xana_va)
117  call xmod%get_field('ud', xmod_ud)
118  call xmod%get_field('vd', xmod_vd)
119  call a2d(geom, xana_ua, xana_va, xmod_ud, xmod_vd)
120  xmod_ud(:,geom%jec+1,:) = 0.0_kind_real
121  xmod_vd(geom%iec+1,:,:) = 0.0_kind_real
122  failed = .false.
123  if (xmod%f_comm%rank() == 0) write(*,"(A)") &
124  "A2M Multiply: analysis increment ua => linearized model ud"
125  endif
126 
127  elseif (xmod%fields(index_mod)%fv3jedi_name == 'vd') then
128 
129  !Special case: A-grid analysis, D-Grid model
130  if (xana%has_field('ua')) then
131  !Already done above
132  failed = .false.
133  if (xmod%f_comm%rank() == 0) write(*,"(A)") &
134  "A2M Multiply: analysis increment va => linearized model ud"
135  endif
136 
137  elseif (xmod%fields(index_mod)%fv3jedi_name == 'delp') then
138 
139  !Special case: ps in analysis, delp in model
140  if (xana%has_field('ps')) then
141  call xana%get_field('ps', xana_ps)
142  call xmod%get_field('delp', xmod_delp)
143  do k = 1,geom%npz
144  xmod_delp(:,:,k) = (geom%bk(k+1)-geom%bk(k))*xana_ps(:,:,1)
145  enddo
146  failed = .false.
147  if (xmod%f_comm%rank() == 0) write(*,"(A)") &
148  "A2M Multiply: analysis increment ps => linearized model ud"
149  endif
150 
151  endif
152 
153  if (failed) call abor1_ftn("fv3jedi_linvarcha_a2m_mod.multiply: found no way of getting "//&
154  trim(xmod%fields(index_mod)%fv3jedi_name)//" from the analysis increment" )
155 
156 enddo
157 
158 ! Copy calendar infomation
159 xmod%calendar_type = xana%calendar_type
160 xmod%date_init = xana%date_init
161 
162 end subroutine multiply
163 
164 ! ------------------------------------------------------------------------------
165 
166 subroutine multiplyadjoint(self,geom,xmod,xana)
167 
168 implicit none
169 type(fv3jedi_linvarcha_a2m), intent(inout) :: self
170 type(fv3jedi_geom), intent(inout) :: geom
171 type(fv3jedi_increment), intent(in) :: xmod
172 type(fv3jedi_increment), intent(inout) :: xana
173 
174 integer :: index_ana, index_mod, index_mod_found
175 integer :: k
176 logical :: failed
177 
178 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ua
179 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_va
180 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ps
181 
182 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_ud
183 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_vd
184 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_delp
185 
186 do index_ana = 1, xana%nf
187 
188  index_mod_found = -1
189  failed = .true.
190 
191  !Check analysis for presence of field
192  do index_mod = 1, xmod%nf
193  if (xana%fields(index_ana)%fv3jedi_name == xmod%fields(index_mod)%fv3jedi_name) then
194  index_mod_found = index_mod
195  exit
196  endif
197  enddo
198 
199  if (index_mod_found >= 0) then
200 
201  !OK, direct copy
202  xana%fields(index_ana)%array = xmod%fields(index_mod_found)%array
203  failed = .false.
204  if (xana%f_comm%rank() == 0) write(*,"(A, A10, A, A10)") &
205  "A2M MultiplyAdjoint: linearized model "//xmod%fields(index_mod_found)%fv3jedi_name(1:10)&
206  //" => analysis increment "//xana%fields(index_ana)%fv3jedi_name(1:10)
207 
208  elseif (xana%fields(index_ana)%fv3jedi_name == 'ua') then
209 
210  !Special case: A-grid analysis, D-Grid model
211  if (xmod%has_field('ud')) then
212  call xana%get_field('ua', xana_ua)
213  call xana%get_field('va', xana_va)
214  call xmod%get_field('ud', xmod_ud)
215  call xmod%get_field('vd', xmod_vd)
216  xmod_ud(:,geom%jec+1,:) = 0.0_kind_real
217  xmod_vd(geom%iec+1,:,:) = 0.0_kind_real
218  call a2d_ad(geom, xana_ua, xana_va, xmod_ud, xmod_vd)
219  failed = .false.
220  if (xana%f_comm%rank() == 0) write(*,"(A)") &
221  "A2M MultiplyAdjoint: linearized model ud => analysis increment ua"
222  endif
223 
224  elseif (xana%fields(index_ana)%fv3jedi_name == 'va') then
225 
226  !Special case: A-grid analysis, D-Grid model
227  if (xmod%has_field('ud')) then
228  !Already done above
229  failed = .false.
230  if (xana%f_comm%rank() == 0) write(*,"(A)") &
231  "A2M MultiplyAdjoint: linearized model vd => analysis increment va"
232  endif
233 
234  elseif (xana%fields(index_ana)%fv3jedi_name == 'ps') then
235 
236  !Special case: ps in analysis, delp in model
237  if (xmod%has_field('delp')) then
238  call xana%get_field('ps', xana_ps)
239  call xmod%get_field('delp', xmod_delp)
240  xana_ps = 0.0_kind_real
241  do k = 1,geom%npz
242  xana_ps(:,:,1) = xana_ps(:,:,1) + (geom%bk(k+1)-geom%bk(k))*xmod_delp(:,:,k)
243  enddo
244  failed = .false.
245  if (xana%f_comm%rank() == 0) write(*,"(A)") &
246  "A2M MultiplyAdjoint: linearized model delp => analysis increment ps"
247  endif
248 
249  endif
250 
251  if (failed) call abor1_ftn("fv3jedi_linvarcha_a2m_mod.multiplyadjoint: found no way of getting "//&
252  trim(xana%fields(index_ana)%fv3jedi_name)//" from the linearized model" )
253 
254 enddo
255 
256 ! Copy calendar infomation
257 xana%calendar_type = xmod%calendar_type
258 xana%date_init = xmod%date_init
259 
260 end subroutine multiplyadjoint
261 
262 ! ------------------------------------------------------------------------------
263 
264 subroutine multiplyinverse(self,geom,xmod,xana)
265 
266 implicit none
267 type(fv3jedi_linvarcha_a2m), intent(inout) :: self
268 type(fv3jedi_geom), intent(inout) :: geom
269 type(fv3jedi_increment), intent(in) :: xmod
270 type(fv3jedi_increment), intent(inout) :: xana
271 
272 integer :: index_ana, index_mod, index_mod_found
273 logical :: failed
274 
275 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ua
276 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_va
277 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ps
278 
279 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_ud
280 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_vd
281 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_delp
282 
283 do index_ana = 1, xana%nf
284 
285  index_mod_found = -1
286  failed = .true.
287 
288  !Check analysis for presence of field
289  do index_mod = 1, xmod%nf
290  if (xana%fields(index_ana)%fv3jedi_name == xmod%fields(index_mod)%fv3jedi_name) then
291  index_mod_found = index_mod
292  exit
293  endif
294  enddo
295 
296  if (index_mod_found >= 0) then
297 
298  !OK, direct copy
299  failed = .false.
300  xana%fields(index_ana)%array = xmod%fields(index_mod_found)%array
301  if (xana%f_comm%rank() == 0) write(*,"(A, A10, A, A10)") &
302  "A2M MultiplyInverse: linearized model "//xmod%fields(index_mod_found)%fv3jedi_name(1:10)&
303  //" => analysis increment "//xana%fields(index_ana)%fv3jedi_name(1:10)
304 
305  elseif (xana%fields(index_ana)%fv3jedi_name == 'ua') then
306 
307  !Special case: A-grid analysis, D-Grid model
308  if (xmod%has_field('ud')) then
309  call xana%get_field('ua', xana_ua)
310  call xana%get_field('va', xana_va)
311  call xmod%get_field('ud', xmod_ud)
312  call xmod%get_field('vd', xmod_vd)
313  xmod_ud(:,geom%jec+1,:) = 0.0_kind_real
314  xmod_vd(geom%iec+1,:,:) = 0.0_kind_real
315  call d2a(geom, xmod_ud, xmod_vd, xana_ua, xana_va)
316  failed = .false.
317  if (xana%f_comm%rank() == 0) write(*,"(A)") &
318  "A2M MultiplyInverse: linearized model ud => analysis increment ua"
319  endif
320 
321  elseif (xana%fields(index_ana)%fv3jedi_name == 'va') then
322 
323  !Special case: A-grid analysis, D-Grid model
324  if (xmod%has_field('ud')) then
325  !Already done above
326  failed = .false.
327  if (xana%f_comm%rank() == 0) write(*,"(A)") &
328  "A2M MultiplyInverse: linearized model vd => analysis increment va"
329  endif
330 
331  elseif (xana%fields(index_ana)%fv3jedi_name == 'ps') then
332 
333  !Special case: ps in analysis, delp in model
334  if (xmod%has_field('delp')) then
335  call xana%get_field('ps', xana_ps)
336  call xmod%get_field('delp', xmod_delp)
337  xana_ps(:,:,1) = sum(xmod_delp,3)
338  failed = .false.
339  if (xana%f_comm%rank() == 0) write(*,"(A)") &
340  "A2M MultiplyInverse: linearized model delp => analysis increment ps"
341  endif
342 
343  endif
344 
345  if (failed) call abor1_ftn("fv3jedi_linvarcha_a2m_mod.multiplyinverse: found no way of getting "//&
346  trim(xana%fields(index_ana)%fv3jedi_name)//" from the linearized model" )
347 
348 enddo
349 
350 ! Copy calendar infomation
351 xana%calendar_type = xmod%calendar_type
352 xana%date_init = xmod%date_init
353 
354 end subroutine multiplyinverse
355 
356 ! ------------------------------------------------------------------------------
357 
358 subroutine multiplyinverseadjoint(self,geom,xana,xmod)
359 
360 implicit none
361 type(fv3jedi_linvarcha_a2m), intent(inout) :: self
362 type(fv3jedi_geom), intent(inout) :: geom
363 type(fv3jedi_increment), intent(in) :: xana
364 type(fv3jedi_increment), intent(inout) :: xmod
365 
366 integer :: index_ana, index_mod, index_ana_found
367 integer :: k
368 logical :: failed
369 
370 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ua
371 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_va
372 real(kind=kind_real), pointer, dimension(:,:,:) :: xana_ps
373 
374 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_ud
375 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_vd
376 real(kind=kind_real), pointer, dimension(:,:,:) :: xmod_delp
377 
378 do index_mod = 1, xmod%nf
379 
380  index_ana_found = -1
381  failed = .true.
382 
383  !Check analysis for presence of field
384  do index_ana = 1, xana%nf
385  if (xmod%fields(index_mod)%fv3jedi_name == xana%fields(index_ana)%fv3jedi_name) then
386  index_ana_found = index_ana
387  exit
388  endif
389  enddo
390 
391  if (index_ana_found >= 0) then
392 
393  !OK, direct copy
394  failed = .false.
395  xmod%fields(index_mod)%array = xana%fields(index_ana_found)%array
396  if (xmod%f_comm%rank() == 0) write(*,"(A, A10, A, A10)") &
397  "A2M MultiplyInverseAdjoint: analysis increment "//xana%fields(index_ana_found)%fv3jedi_name(1:10)&
398  //" => linearized model "//xmod%fields(index_mod)%fv3jedi_name(1:10)
399 
400  elseif (xmod%fields(index_mod)%fv3jedi_name == 'ud') then
401 
402  !Special case: A-grid analysis, D-Grid model
403  if (xana%has_field('ua')) then
404  call xana%get_field('ua', xana_ua)
405  call xana%get_field('va', xana_va)
406  call xmod%get_field('ud', xmod_ud)
407  call xmod%get_field('vd', xmod_vd)
408  call d2a_ad(geom, xmod_ud, xmod_vd, xana_ua, xana_va)
409  xmod_ud(:,geom%jec+1,:) = 0.0_kind_real
410  xmod_vd(geom%iec+1,:,:) = 0.0_kind_real
411  failed = .false.
412  if (xmod%f_comm%rank() == 0) write(*,"(A)") &
413  "A2M MultiplyInverseAdjoint: analysis increment ua => linearized model ud"
414  endif
415 
416  elseif (xmod%fields(index_mod)%fv3jedi_name == 'vd') then
417 
418  !Special case: A-grid analysis, D-Grid model
419  if (xana%has_field('ua')) then
420  !Already done above
421  failed = .false.
422  if (xmod%f_comm%rank() == 0) write(*,"(A)") &
423  "A2M MultiplyInverseAdjoint: analysis increment va => linearized model vd"
424  endif
425 
426  elseif (xmod%fields(index_mod)%fv3jedi_name == 'delp') then
427 
428  !Special case: ps in analysis, delp in model
429  if (xana%has_field('ps')) then
430  call xana%get_field('ps', xana_ps)
431  call xmod%get_field('delp', xmod_delp)
432  do k = 1,geom%npz
433  xmod_delp(:,:,k) = xana_ps(:,:,1)
434  enddo
435  failed = .false.
436  if (xmod%f_comm%rank() == 0) write(*,"(A)") &
437  "A2M MultiplyInverseAdjoint: analysis increment ps => linearized model delp"
438  endif
439 
440  endif
441 
442  if (failed) call abor1_ftn("fv3jedi_linvarcha_a2m_mod.multiplyinverseadjoint: found no way of getting "//&
443  trim(xmod%fields(index_mod)%fv3jedi_name)//" from the analysis increment" )
444 
445 enddo
446 
447 ! Copy calendar infomation
448 xmod%calendar_type = xana%calendar_type
449 xmod%date_init = xana%date_init
450 
451 end subroutine multiplyinverseadjoint
452 
453 ! ------------------------------------------------------------------------------
454 
455 end module fv3jedi_linvarcha_a2m_mod
fv3jedi_state_mod::fv3jedi_state
Fortran derived type to hold FV3JEDI state.
Definition: fv3jedi_state_mod.F90:30
fv3jedi_field_mod
Definition: fv3jedi_field_mod.f90:6
fv3jedi_linvarcha_a2m_mod::multiply
subroutine, public multiply(self, geom, xana, xmod)
Definition: fv3jedi_linvarcha_a2m_mod.f90:70
fv3jedi_linvarcha_a2m_mod::multiplyinverseadjoint
subroutine, public multiplyinverseadjoint(self, geom, xana, xmod)
Definition: fv3jedi_linvarcha_a2m_mod.f90:359
fv3jedi_state_mod
Definition: fv3jedi_state_mod.F90:6
fv3jedi_field_mod::copy_subset
subroutine, public copy_subset(field_in, field_ou, not_copied)
Definition: fv3jedi_field_mod.f90:236
fv3jedi_geom_mod
Fortran module handling geometry for the FV3 model.
Definition: fv3jedi_geom_mod.f90:8
fv3jedi_increment_mod
Definition: fv3jedi_increment_mod.F90:6
fv3jedi_geom_mod::fv3jedi_geom
Fortran derived type to hold geometry data for the FV3JEDI model.
Definition: fv3jedi_geom_mod.f90:46
wind_vt_mod::d2a_ad
subroutine, public d2a_ad(geom, u_ad_comp, v_ad_comp, ua_ad_comp, va_ad_comp)
Definition: wind_variables_mod.f90:1607
wind_vt_mod
Definition: wind_variables_mod.f90:6
fv3jedi_linvarcha_a2m_mod
Definition: fv3jedi_linvarcha_a2m_mod.f90:6
fv3jedi_linvarcha_a2m_mod::multiplyadjoint
subroutine, public multiplyadjoint(self, geom, xmod, xana)
Definition: fv3jedi_linvarcha_a2m_mod.f90:167
fv3jedi_linvarcha_a2m_mod::delete
subroutine, public delete(self)
Definition: fv3jedi_linvarcha_a2m_mod.f90:61
fv3jedi_linvarcha_a2m_mod::fv3jedi_linvarcha_a2m
Definition: fv3jedi_linvarcha_a2m_mod.f90:31
fv3jedi_kinds_mod::kind_real
integer, parameter, public kind_real
Definition: fv3jedi_kinds_mod.f90:14
wind_vt_mod::d2a
subroutine, public d2a(geom, u_comp, v_comp, ua_comp, va_comp)
Definition: wind_variables_mod.f90:1456
wind_vt_mod::a2d
subroutine, public a2d(geom, ua, va, ud, vd)
Definition: wind_variables_mod.f90:1023
fv3jedi_linvarcha_a2m_mod::multiplyinverse
subroutine, public multiplyinverse(self, geom, xmod, xana)
Definition: fv3jedi_linvarcha_a2m_mod.f90:265
fv3jedi_kinds_mod
Definition: fv3jedi_kinds_mod.f90:6
wind_vt_mod::a2d_ad
subroutine, public a2d_ad(geom, ua_ad, va_ad, ud_ad, vd_ad)
Definition: wind_variables_mod.f90:1194
fv3jedi_increment_mod::fv3jedi_increment
Definition: fv3jedi_increment_mod.F90:34
fv3jedi_linvarcha_a2m_mod::create
subroutine, public create(self, geom, bg, fg, c_conf)
Definition: fv3jedi_linvarcha_a2m_mod.f90:42