Interpolate 3D FV data into DG nodes.
62 use scale_comm_cartesc, only: &
63 comm_vars8, &
64 comm_wait
65 implicit none
66 class(LocalMesh3D), intent(in) :: lcmesh3D
67 class(ElementBase3D), intent(in) :: elem3D_
68 integer, intent(in) :: IS, IE, IA, IHALO
69 integer, intent(in) :: JS, JE, JA, JHALO
70 integer, intent(in) :: KS, KE, KA, KHALO
71 real(RP), intent(inout) :: out_dg(elem3D_%Np,lcmesh3D%NeA)
72 real(RP), intent(in) :: in_fv(KA,IA,JA)
73 integer, intent(in) :: interp_ord
74 logical, intent(in), optional :: out_dg_flush_flag
75
76 integer :: ke_x, ke_y, ke_z
77 integer :: i, j, k
78 integer :: kelem
79
80 real(RP) :: tend_fv_cp(KA,IA,JA)
81 real(RP) :: dqdx(KA,IA,JA), dqdy(KA,IA,JA), dqdz(KA,IA,JA)
82 real(RP) :: dqdxx(KA,IA,JA), dqdyy(KA,IA,JA), dqdzz(KA,IA,JA)
83 real(RP) :: dqdxy(KA,IA,JA), dqdxz(KA,IA,JA), dqdyz(KA,IA,JA)
84
85 logical :: out_dg_flush_flag_
86
87
88 if ( present(out_dg_flush_flag) ) then
89 out_dg_flush_flag_ = out_dg_flush_flag
90 else
91 out_dg_flush_flag_ = .true.
92 end if
93
94 tend_fv_cp(:,:,:) = in_fv(:,:,:)
95 if ( interp_ord > 0 ) then
96 call comm_vars8( tend_fv_cp(:,:,:), 1 )
97 call comm_wait( tend_fv_cp(:,:,:), 1, .true. )
98 end if
99
100
101
102 if ( out_dg_flush_flag_ ) then
103
104 out_dg(:,:) = 0.0_rp
105
106 end if
107
108 if ( interp_ord > 0 ) then
109
110 do j=js-1, je+1
111 do i=is, ie
112 do k=ks, ke
113 dqdx(k,i,j) = 0.25_rp * ( tend_fv_cp(k,i+1,j) - tend_fv_cp(k,i-1,j) )
114 enddo
115 enddo
116 enddo
117
118 do j=js, je
119 do i=is-1, ie+1
120 do k=ks, ke
121 dqdy(k,i,j) = 0.25_rp * ( tend_fv_cp(k,i,j+1) - tend_fv_cp(k,i,j-1) )
122 enddo
123 enddo
124 enddo
125
126 do j=js, je
127 do i=is, ie
128 do k=ks+1, ke-1
129 dqdz(k,i,j) = 0.25_rp * ( tend_fv_cp(k+1,i,j) - tend_fv_cp(k-1,i,j) )
130 enddo
131 dqdz(ks,i,j) = 0.25_rp * ( - tend_fv_cp(ks+2,i,j) + 4.0_rp * tend_fv_cp(ks+1,i,j) - 3.0_rp * tend_fv_cp(ks,i,j) )
132 dqdz(ke,i,j) = 0.25_rp * ( 3.0_rp * tend_fv_cp(ke,i,j) - 4.0_rp * tend_fv_cp(ke-1,i,j) + tend_fv_cp(ke-2,i,j) )
133 enddo
134 enddo
135 end if
136 if ( interp_ord > 1 ) then
137
138 do j=js, je
139 do i=is, ie
140 do k=ks, ke
141 dqdxx(k,i,j) = 0.25_rp * ( tend_fv_cp(k,i+1,j) - 2.0_rp * tend_fv_cp(k,i,j) + tend_fv_cp(k,i-1,j) )
142 enddo
143 enddo
144 enddo
145
146 do j=js, je
147 do i=is, ie
148 do k=ks, ke
149 dqdyy(k,i,j) = 0.25_rp * ( tend_fv_cp(k,i,j+1) - 2.0_rp * tend_fv_cp(k,i,j) + tend_fv_cp(k,i,j-1) )
150 enddo
151 enddo
152 enddo
153
154 do j=js, je
155 do i=is, ie
156 do k=ks, ke
157 dqdxy(k,i,j) = 0.25_rp * ( dqdx(k,i,j+1) - dqdx(k,i,j-1) )
158 enddo
159 enddo
160 enddo
161
162 do j=js, je
163 do i=is, ie
164 do k=ks+1, ke-1
165 dqdzz(k,i,j) = 0.25_rp * ( tend_fv_cp(k+1,i,j) - 2.0_rp * tend_fv_cp(k,i,j) + tend_fv_cp(k-1,i,j) )
166 dqdxz(k,i,j) = 0.25_rp * ( dqdx(k+1,i,j) - dqdx(k-1,i,j) )
167 dqdyz(k,i,j) = 0.25_rp * ( dqdy(k+1,i,j) - dqdy(k-1,i,j) )
168 enddo
169 dqdzz(ks,i,j) = dqdzz(ks+1,i,j); dqdzz(ke,i,j) = dqdzz(ke-1,i,j);
170 dqdxz(ks,i,j) = 0.5_rp * ( dqdx(ks+1,i,j) - dqdx(ks,i,j) ); dqdxz(ke,i,j) = 0.5_rp * ( dqdx(ke,i,j) - dqdx(ke-1,i,j) );
171 dqdyz(ks,i,j) = 0.5_rp * ( dqdy(ks+1,i,j) - dqdy(ks,i,j) ); dqdyz(ke,i,j) = 0.5_rp * ( dqdy(ke,i,j) - dqdy(ke-1,i,j) );
172 enddo
173 enddo
174 end if
175
176 select case(interp_ord)
177 case(0)
178
179 do ke_z=1, lcmesh3d%NeZ
180 do ke_y=1, lcmesh3d%NeY
181 do ke_x=1, lcmesh3d%NeX
182 kelem = ke_x + (ke_y-1)*lcmesh3d%NeX + (ke_z-1)*lcmesh3d%NeX*lcmesh3d%NeY
183 i = ihalo + ke_x; j = jhalo + ke_y; k = khalo + ke_z
184 out_dg(:,kelem) = out_dg(:,kelem) &
185 + in_fv(k,i,j)
186 enddo
187 enddo
188 enddo
189 case(1)
190
191 do ke_z=1, lcmesh3d%NeZ
192 do ke_y=1, lcmesh3d%NeY
193 do ke_x=1, lcmesh3d%NeX
194 kelem = ke_x + (ke_y-1)*lcmesh3d%NeX + (ke_z-1)*lcmesh3d%NeX*lcmesh3d%NeY
195 i = ihalo + ke_x; j = jhalo + ke_y; k = khalo + ke_z
196 out_dg(:,kelem) = out_dg(:,kelem) &
197 + in_fv(k,i,j) &
198 + dqdx(k,i,j) * elem3d_%x1(:) &
199 + dqdy(k,i,j) * elem3d_%x2(:) &
200 + dqdz(k,i,j) * elem3d_%x3(:)
201 enddo
202 enddo
203 enddo
204 case(2)
205
206 do ke_z=1, lcmesh3d%NeZ
207 do ke_y=1, lcmesh3d%NeY
208 do ke_x=1, lcmesh3d%NeX
209 kelem = ke_x + (ke_y-1)*lcmesh3d%NeX + (ke_z-1)*lcmesh3d%NeX*lcmesh3d%NeY
210 i = ihalo + ke_x; j = jhalo + ke_y; k = khalo + ke_z
211 out_dg(:,kelem) = out_dg(:,kelem) &
212 + in_fv(k,i,j) &
213 + dqdx(k,i,j) * elem3d_%x1(:) &
214 + dqdy(k,i,j) * elem3d_%x2(:) &
215 + dqdz(k,i,j) * elem3d_%x3(:) &
216 + dqdxy(k,i,j) * elem3d_%x1(:) * elem3d_%x2(:) &
217 + dqdxz(k,i,j) * elem3d_%x1(:) * elem3d_%x3(:) &
218 + dqdyz(k,i,j) * elem3d_%x2(:) * elem3d_%x3(:) &
219 + 0.5_rp * dqdxx(k,i,j) * ( elem3d_%x1(:)**2 - 1.0_rp / 3.0_rp ) &
220 + 0.5_rp * dqdyy(k,i,j) * ( elem3d_%x2(:)**2 - 1.0_rp / 3.0_rp ) &
221 + 0.5_rp * dqdzz(k,i,j) * ( elem3d_%x3(:)**2 - 1.0_rp / 3.0_rp )
222 enddo
223 enddo
224 enddo
225 end select
226
227 return