CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
mp_thompson_post.F90
1
3
6
7 use mpi_f08
8 use machine, only : kind_phys
9
10 implicit none
11
12 public :: mp_thompson_post_init, mp_thompson_post_run, mp_thompson_post_finalize
13
14 private
15
16 logical :: is_initialized = .false.
17
18 logical :: apply_limiter
19
20contains
21
25 subroutine mp_thompson_post_init(ttendlim, errmsg, errflg)
26
27 implicit none
28
29 ! Interface variables
30 real(kind_phys), intent(in) :: ttendlim
31
32 ! CCPP error handling
33 character(len=*), intent( out) :: errmsg
34 integer, intent( out) :: errflg
35
36 ! Local variables
37 integer :: i
38
39 ! Initialize the CCPP error handling variables
40 errmsg = ''
41 errflg = 0
42
43 ! Check initialization state
44 if (is_initialized) return
45
46 if (ttendlim < 0) then
47 apply_limiter = .false.
48 else
49 apply_limiter = .true.
50 end if
51
52 is_initialized = .true.
53
54 end subroutine mp_thompson_post_init
55
59 subroutine mp_thompson_post_run(ncol, nlev, tgrs_save, tgrs, prslk, dtp, ttendlim, &
60 kdt, mpicomm, mpirank, mpiroot, errmsg, errflg)
61
62 implicit none
63
64 ! Interface variables
65 integer, intent(in) :: ncol
66 integer, intent(in) :: nlev
67 real(kind_phys), dimension(:,:), intent(in) :: tgrs_save
68 real(kind_phys), dimension(:,:), intent(inout) :: tgrs
69 real(kind_phys), dimension(:,:), intent(in) :: prslk
70 real(kind_phys), intent(in) :: dtp
71 real(kind_phys), intent(in) :: ttendlim
72 integer, intent(in) :: kdt
73 ! MPI information
74 type(mpi_comm), intent(in ) :: mpicomm
75 integer, intent(in ) :: mpirank
76 integer, intent(in ) :: mpiroot
77 ! CCPP error handling
78 character(len=*), intent( out) :: errmsg
79 integer, intent( out) :: errflg
80
81 ! Local variables
82 real(kind_phys), dimension(1:ncol,1:nlev) :: mp_tend
83 integer :: i, k
84#ifdef DEBUG
85 integer :: events
86#endif
87
88 ! Initialize the CCPP error handling variables
89 errmsg = ''
90 errflg = 0
91
92 ! Check initialization state
93 if (.not.is_initialized) then
94 write(errmsg, fmt='((a))') 'mp_thompson_post_run called before mp_thompson_post_init'
95 errflg = 1
96 return
97 end if
98
99 ! If limiter is deactivated, return immediately
100 if (.not.apply_limiter) return
101
102 ! mp_tend and ttendlim are expressed in potential temperature
103 mp_tend = (tgrs - tgrs_save)/prslk
104
105#ifdef DEBUG
106 events = 0
107#endif
108 do k=1,nlev
109 do i=1,ncol
110 mp_tend(i,k) = max( -ttendlim*dtp, min( ttendlim*dtp, mp_tend(i,k) ) )
111
112#ifdef DEBUG
113 if (tgrs_save(i,k) + mp_tend(i,k)*prslk(i,k) .ne. tgrs(i,k)) then
114 write(0,'(a,3i6,3e16.7)') "mp_thompson_post_run mp_tend limiter: kdt, i, k, t_old, t_new, t_lim:", &
115 & kdt, i, k, tgrs_save(i,k), tgrs(i,k), tgrs_save(i,k) + mp_tend(i,k)*prslk(i,k)
116 events = events + 1
117 end if
118#endif
119 tgrs(i,k) = tgrs_save(i,k) + mp_tend(i,k)*prslk(i,k)
120 end do
121 end do
122
123#ifdef DEBUG
124 if (events > 0) then
125 write(0,'(a,i0,a,i0,a,i0)') "mp_thompson_post_run: ttendlim applied ", events, "/", nlev*ncol, &
126 & " times at timestep ", kdt
127 end if
128#endif
129
130 end subroutine mp_thompson_post_run
131
135 subroutine mp_thompson_post_finalize(errmsg, errflg)
136
137 implicit none
138
139 ! CCPP error handling
140 character(len=*), intent( out) :: errmsg
141 integer, intent( out) :: errflg
142
143 ! initialize ccpp error handling variables
144 errmsg = ''
145 errflg = 0
146
147 ! Check initialization state
148 if (.not. is_initialized) return
149
150 is_initialized = .false.
151
152 end subroutine mp_thompson_post_finalize
153
154end module mp_thompson_post
This module contain the post processing of Thompson microphysics.