FV3 Bundle
qg_trajectories.F90
Go to the documentation of this file.
1 ! (C) Copyright 2009-2016 ECMWF.
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 ! In applying this licence, ECMWF does not waive the privileges and immunities
6 ! granted to it by virtue of its status as an intergovernmental organisation nor
7 ! does it submit to any jurisdiction.
8 
9 !> Handle the trajectory for the QG model
10 
12 
13 use kinds
14 implicit none
15 private
16 
18 public :: qg_traj_registry
19 
20 ! ------------------------------------------------------------------------------
21 
22 !> Fortran derived type to hold the QG model trajectory
24  real(kind=kind_real), allocatable :: x(:,:,:)
25  real(kind=kind_real), allocatable :: xn(:)
26  real(kind=kind_real), allocatable :: xs(:)
27  real(kind=kind_real), allocatable :: qn(:,:)
28  real(kind=kind_real), allocatable :: qs(:,:)
29 end type qg_trajectory
30 
31 #define LISTED_TYPE qg_trajectory
32 
33 !> Linked list interface - defines registry_t type
34 #include "oops/util/linkedList_i.f"
35 
36 !> Global registry
37 type(registry_t) :: qg_traj_registry
38 
39 ! ------------------------------------------------------------------------------
40 contains
41 ! ------------------------------------------------------------------------------
42 !> Linked list implementation
43 #include "oops/util/linkedList_c.f"
44 
45 ! ------------------------------------------------------------------------------
46 
47 subroutine set_traj(self,kx,ky,px,pxn,pxs,pqn,pqs)
48 implicit none
49 type(qg_trajectory), intent(inout) :: self
50 integer, intent(in) :: kx,ky
51 real(kind=kind_real), intent(in) :: px(kx,ky,2)
52 real(kind=kind_real), intent(in) :: pxn(2),pxs(2)
53 real(kind=kind_real), intent(in) :: pqn(kx,2),pqs(kx,2)
54 
55 allocate(self%x(kx,ky,2))
56 allocate(self%xn(2),self%xs(2))
57 allocate(self%qn(kx,2),self%qs(kx,2))
58 
59 self%x(:,:,:)=px(:,:,:)
60 self%xn(:) =pxn(:)
61 self%xs(:) =pxs(:)
62 self%qn(:,:) =pqn(:,:)
63 self%qs(:,:) =pqs(:,:)
64 
65 return
66 end subroutine set_traj
67 
68 ! ------------------------------------------------------------------------------
69 
70 subroutine get_traj(self,kx,ky,px,pxn,pxs,pqn,pqs)
71 implicit none
72 type(qg_trajectory), intent(in) :: self
73 integer, intent(in) :: kx,ky
74 real(kind=kind_real), intent(inout) :: px(kx,ky,2)
75 real(kind=kind_real), intent(inout) :: pxn(2),pxs(2)
76 real(kind=kind_real), intent(inout) :: pqn(kx,2),pqs(kx,2)
77 
78 px(:,:,:)=self%x(:,:,:)
79 pxn(:) =self%xn(:)
80 pxs(:) =self%xs(:)
81 pqn(:,:) =self%qn(:,:)
82 pqs(:,:) =self%qs(:,:)
83 
84 return
85 end subroutine get_traj
86 
87 ! ------------------------------------------------------------------------------
88 
89 subroutine delete_traj(self)
90 implicit none
91 type(qg_trajectory), intent(inout) :: self
92 
93 deallocate(self%x)
94 deallocate(self%xn)
95 deallocate(self%xs)
96 deallocate(self%qn)
97 deallocate(self%qs)
98 
99 return
100 end subroutine delete_traj
101 
102 ! ------------------------------------------------------------------------------
103 
104 subroutine c_minmax_traj(c_key_self, pminmax) bind(c,name='qg_traj_minmaxrms_f90')
105 use iso_c_binding
106 implicit none
107 integer(c_int), intent(in) :: c_key_self
108 real(c_double), intent(inout) :: pminmax(3,5)
109 type(qg_trajectory), pointer :: self
110 real(kind=kind_real) :: zz
111 
112 call qg_traj_registry%get(c_key_self,self)
113 
114 zz=real(size(self%x),kind_real)
115 pminmax(1,1)=minval(self%x(:,:,:))
116 pminmax(2,1)=maxval(self%x(:,:,:))
117 pminmax(3,1)=sqrt(sum(self%x(:,:,:)**2)/zz)
118 
119 zz=2.0_kind_real
120 pminmax(1,2)=minval(self%xn(:))
121 pminmax(2,2)=maxval(self%xn(:))
122 pminmax(3,2)=sqrt(sum(self%xn(:)**2)/zz)
123 
124 pminmax(1,3)=minval(self%xs(:))
125 pminmax(2,3)=maxval(self%xs(:))
126 pminmax(3,3)=sqrt(sum(self%xs(:)**2)/zz)
127 
128 zz=real(size(self%qn),kind_real)
129 pminmax(1,4)=minval(self%qn(:,:))
130 pminmax(2,4)=maxval(self%qn(:,:))
131 pminmax(3,4)=sqrt(sum(self%qn(:,:)**2)/zz)
132 
133 pminmax(1,5)=minval(self%qs(:,:))
134 pminmax(2,5)=maxval(self%qs(:,:))
135 pminmax(3,5)=sqrt(sum(self%qs(:,:)**2)/zz)
136 
137 end subroutine c_minmax_traj
138 
139 ! ------------------------------------------------------------------------------
140 
141 end module qg_trajectories
Handle the trajectory for the QG model.
subroutine, public set_traj(self, kx, ky, px, pxn, pxs, pqn, pqs)
Linked list implementation.
subroutine, public get_traj(self, kx, ky, px, pxn, pxs, pqn, pqs)
type(registry_t), public qg_traj_registry
Linked list interface - defines registry_t type.
Fortran derived type to hold the QG model trajectory.
subroutine c_minmax_traj(c_key_self, pminmax)
subroutine, public delete_traj(self)