FV3 Bundle
fv3jedi_traj_mod.f90
Go to the documentation of this file.
1 ! (C) Copyright 2017 UCAR
2 !
3 ! This software is licensed under the terms of the Apache Licence Version 2.0
4 ! which can be obtained at http://www.apache.org/licenses/LICENSE-2.0.
5 
7 
9 use iso_c_binding
12 use fv3jedi_lm_utils_mod, only: fv3jedi_traj => fv3jedi_lm_traj, &
14 
15 implicit none
16 private
17 
18 public :: fv3jedi_traj
19 public :: traj_prop
20 public :: traj_wipe
21 public :: traj_minmaxrms
22 
23 ! ------------------------------------------------------------------------------
24 
25 contains
26 
27 ! ------------------------------------------------------------------------------
28 
29 subroutine traj_prop(model, state, self)
30 
31 implicit none
32 type(fv3jedi_model) :: model
33 type(fv3jedi_state) :: state
34 type(fv3jedi_traj) :: self
35 
36 call allocate_traj(self,state%isc,state%iec,state%jsc,state%jec,state%npz,&
37  model%fv3jedi_lm%conf%hydrostatic,model%fv3jedi_lm%conf%do_phy_mst)
38 
39 self%u = state%ud
40 self%v = state%vd
41 self%ua = state%ua
42 self%va = state%va
43 self%t = state%t
44 self%delp = state%delp
45 self%qv = state%q
46 self%ql = state%qi
47 self%qi = state%ql
48 self%o3 = state%o3
49 
50 if (.not. model%fv3jedi_lm%conf%hydrostatic) then
51 self%w = state%w
52 self%delz = state%delz
53 endif
54 
55 if (model%fv3jedi_lm%conf%do_phy_mst /= 0) then
56 self%qls = state%qls
57 self%qcn = state%qcn
58 self%cfcn = state%cfcn
59 endif
60 
61 self%phis = state%phis
62 self%frocean = state%frocean
63 self%frland = state%frland
64 self%varflt = state%varflt
65 self%ustar = state%ustar
66 self%bstar = state%bstar
67 self%zpbl = state%zpbl
68 self%cm = state%cm
69 self%ct = state%ct
70 self%cq = state%cq
71 self%kcbl = state%kcbl
72 self%ts = state%ts
73 self%khl = state%khl
74 self%khu = state%khu
75 
76 end subroutine traj_prop
77 
78 ! ------------------------------------------------------------------------------
79 
80 subroutine traj_wipe(self)
81 
82 implicit none
83 type(fv3jedi_traj), pointer :: self
84 
85 call deallocate_traj(self)
86 
87 end subroutine traj_wipe
88 
89 ! ------------------------------------------------------------------------------
90 
91 subroutine traj_minmaxrms(self, pminmax)
92 
93 implicit none
94 type(fv3jedi_traj), intent(in) :: self
95 real(c_double), intent(inout) :: pminmax(3,11)
96 
97 real(kind=kind_real) :: zz
98 
99 pminmax = 0.0_kind_real
100 
101 zz=real(size(self%u),kind_real)
102 pminmax(1,1)=minval(self%u(:,:,:))
103 pminmax(2,1)=maxval(self%u(:,:,:))
104 pminmax(3,1)=sqrt(sum(self%u(:,:,:)**2)/zz)
105 
106 zz=real(size(self%v),kind_real)
107 pminmax(1,2)=minval(self%v(:,:,:))
108 pminmax(2,2)=maxval(self%v(:,:,:))
109 pminmax(3,2)=sqrt(sum(self%v(:,:,:)**2)/zz)
110 
111 zz=real(size(self%t),kind_real)
112 pminmax(1,3)=minval(self%t(:,:,:))
113 pminmax(2,3)=maxval(self%t(:,:,:))
114 pminmax(3,3)=sqrt(sum(self%t(:,:,:)**2)/zz)
115 
116 zz=real(size(self%delp),kind_real)
117 pminmax(1,4)=minval(self%delp(:,:,:))
118 pminmax(2,4)=maxval(self%delp(:,:,:))
119 pminmax(3,4)=sqrt(sum(self%delp(:,:,:)**2)/zz)
120 
121 zz=real(size(self%qv),kind_real)
122 pminmax(1,5)=minval(self%qv(:,:,:))
123 pminmax(2,5)=maxval(self%qv(:,:,:))
124 pminmax(3,5)=sqrt(sum(self%qv(:,:,:)**2)/zz)
125 
126 zz=real(size(self%qi),kind_real)
127 pminmax(1,6)=minval(self%qi(:,:,:))
128 pminmax(2,6)=maxval(self%qi(:,:,:))
129 pminmax(3,6)=sqrt(sum(self%qi(:,:,:)**2)/zz)
130 
131 zz=real(size(self%ql),kind_real)
132 pminmax(1,7)=minval(self%ql(:,:,:))
133 pminmax(2,7)=maxval(self%ql(:,:,:))
134 pminmax(3,7)=sqrt(sum(self%ql(:,:,:)**2)/zz)
135 
136 zz=real(size(self%o3),kind_real)
137 pminmax(1,8)=minval(self%o3(:,:,:))
138 pminmax(2,8)=maxval(self%o3(:,:,:))
139 pminmax(3,8)=sqrt(sum(self%o3(:,:,:)**2)/zz)
140 
141 if (allocated(self%w)) then
142  zz=real(size(self%w),kind_real)
143  pminmax(1,9)=minval(self%w(:,:,:))
144  pminmax(2,9)=maxval(self%w(:,:,:))
145  pminmax(3,9)=sqrt(sum(self%w(:,:,:)**2)/zz)
146 endif
147 
148 if (allocated(self%delz)) then
149  zz=real(size(self%delz),kind_real)
150  pminmax(1,10)=minval(self%delz(:,:,:))
151  pminmax(2,10)=maxval(self%delz(:,:,:))
152  pminmax(3,10)=sqrt(sum(self%delz(:,:,:)**2)/zz)
153 endif
154 
155 zz=real(size(self%phis),kind_real)
156 pminmax(1,11)=minval(self%phis(:,:))
157 pminmax(2,11)=maxval(self%phis(:,:))
158 pminmax(3,11)=sqrt(sum(self%phis(:,:)**2)/zz)
159 
160 end subroutine traj_minmaxrms
161 
162 ! ------------------------------------------------------------------------------
163 
164 end module fv3jedi_traj_mod
subroutine, public allocate_traj(traj, isc, iec, jsc, jec, npz, hydrostatic, dpm)
Fortran derived type to hold FV3JEDI state.
subroutine, public deallocate_traj(traj)
subroutine, public traj_wipe(self)
subroutine, public traj_prop(model, state, self)
subroutine, public traj_minmaxrms(self, pminmax)
Handle state for the FV3JEDI odel.
Fortran derived type to hold the linearized model trajectory.
Fortran derived type to hold model definition.
integer, parameter, public kind_real