9 use fckit_configuration_module,
only: fckit_configuration
41 subroutine create(self, geom, bg, fg, c_conf)
48 type(c_ptr),
intent(in) :: c_conf
50 type(fckit_configuration) :: f_conf
52 f_conf = fckit_configuration(c_conf)
77 integer :: index_ana, index_mod, index_ana_found
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
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
89 do index_mod = 1, xmod%nf
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
102 if (index_ana_found >= 0)
then
105 xmod%fields(index_mod)%array = xana%fields(index_ana_found)%array
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)
111 elseif (xmod%fields(index_mod)%fv3jedi_name ==
'ud')
then
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
123 if (xmod%f_comm%rank() == 0)
write(*,
"(A)") &
124 "A2M Multiply: analysis increment ua => linearized model ud"
127 elseif (xmod%fields(index_mod)%fv3jedi_name ==
'vd')
then
130 if (xana%has_field(
'ua'))
then
133 if (xmod%f_comm%rank() == 0)
write(*,
"(A)") &
134 "A2M Multiply: analysis increment va => linearized model ud"
137 elseif (xmod%fields(index_mod)%fv3jedi_name ==
'delp')
then
140 if (xana%has_field(
'ps'))
then
141 call xana%get_field(
'ps', xana_ps)
142 call xmod%get_field(
'delp', xmod_delp)
144 xmod_delp(:,:,k) = (geom%bk(k+1)-geom%bk(k))*xana_ps(:,:,1)
147 if (xmod%f_comm%rank() == 0)
write(*,
"(A)") &
148 "A2M Multiply: analysis increment ps => linearized model ud"
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" )
159 xmod%calendar_type = xana%calendar_type
160 xmod%date_init = xana%date_init
174 integer :: index_ana, index_mod, index_mod_found
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
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
186 do index_ana = 1, xana%nf
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
199 if (index_mod_found >= 0)
then
202 xana%fields(index_ana)%array = xmod%fields(index_mod_found)%array
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)
208 elseif (xana%fields(index_ana)%fv3jedi_name ==
'ua')
then
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)
220 if (xana%f_comm%rank() == 0)
write(*,
"(A)") &
221 "A2M MultiplyAdjoint: linearized model ud => analysis increment ua"
224 elseif (xana%fields(index_ana)%fv3jedi_name ==
'va')
then
227 if (xmod%has_field(
'ud'))
then
230 if (xana%f_comm%rank() == 0)
write(*,
"(A)") &
231 "A2M MultiplyAdjoint: linearized model vd => analysis increment va"
234 elseif (xana%fields(index_ana)%fv3jedi_name ==
'ps')
then
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
242 xana_ps(:,:,1) = xana_ps(:,:,1) + (geom%bk(k+1)-geom%bk(k))*xmod_delp(:,:,k)
245 if (xana%f_comm%rank() == 0)
write(*,
"(A)") &
246 "A2M MultiplyAdjoint: linearized model delp => analysis increment ps"
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" )
257 xana%calendar_type = xmod%calendar_type
258 xana%date_init = xmod%date_init
272 integer :: index_ana, index_mod, index_mod_found
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
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
283 do index_ana = 1, xana%nf
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
296 if (index_mod_found >= 0)
then
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)
305 elseif (xana%fields(index_ana)%fv3jedi_name ==
'ua')
then
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)
317 if (xana%f_comm%rank() == 0)
write(*,
"(A)") &
318 "A2M MultiplyInverse: linearized model ud => analysis increment ua"
321 elseif (xana%fields(index_ana)%fv3jedi_name ==
'va')
then
324 if (xmod%has_field(
'ud'))
then
327 if (xana%f_comm%rank() == 0)
write(*,
"(A)") &
328 "A2M MultiplyInverse: linearized model vd => analysis increment va"
331 elseif (xana%fields(index_ana)%fv3jedi_name ==
'ps')
then
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)
339 if (xana%f_comm%rank() == 0)
write(*,
"(A)") &
340 "A2M MultiplyInverse: linearized model delp => analysis increment ps"
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" )
351 xana%calendar_type = xmod%calendar_type
352 xana%date_init = xmod%date_init
366 integer :: index_ana, index_mod, index_ana_found
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
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
378 do index_mod = 1, xmod%nf
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
391 if (index_ana_found >= 0)
then
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)
400 elseif (xmod%fields(index_mod)%fv3jedi_name ==
'ud')
then
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
412 if (xmod%f_comm%rank() == 0)
write(*,
"(A)") &
413 "A2M MultiplyInverseAdjoint: analysis increment ua => linearized model ud"
416 elseif (xmod%fields(index_mod)%fv3jedi_name ==
'vd')
then
419 if (xana%has_field(
'ua'))
then
422 if (xmod%f_comm%rank() == 0)
write(*,
"(A)") &
423 "A2M MultiplyInverseAdjoint: analysis increment va => linearized model vd"
426 elseif (xmod%fields(index_mod)%fv3jedi_name ==
'delp')
then
429 if (xana%has_field(
'ps'))
then
430 call xana%get_field(
'ps', xana_ps)
431 call xmod%get_field(
'delp', xmod_delp)
433 xmod_delp(:,:,k) = xana_ps(:,:,1)
436 if (xmod%f_comm%rank() == 0)
write(*,
"(A)") &
437 "A2M MultiplyInverseAdjoint: analysis increment ps => linearized model delp"
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" )
448 xmod%calendar_type = xana%calendar_type
449 xmod%date_init = xana%date_init