Neko  0.9.99
A portable framework for high-order spectral element flow simulations
fluid_pnpn.f90
Go to the documentation of this file.
1 ! Copyright (c) 2022-2024, The Neko Authors
2 ! All rights reserved.
3 !
4 ! Redistribution and use in source and binary forms, with or without
5 ! modification, are permitted provided that the following conditions
6 ! are met:
7 !
8 ! * Redistributions of source code must retain the above copyright
9 ! notice, this list of conditions and the following disclaimer.
10 !
11 ! * Redistributions in binary form must reproduce the above
12 ! copyright notice, this list of conditions and the following
13 ! disclaimer in the documentation and/or other materials provided
14 ! with the distribution.
15 !
16 ! * Neither the name of the authors nor the names of its
17 ! contributors may be used to endorse or promote products derived
18 ! from this software without specific prior written permission.
19 !
20 ! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21 ! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22 ! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23 ! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24 ! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25 ! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26 ! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27 ! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28 ! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29 ! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30 ! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31 ! POSSIBILITY OF SUCH DAMAGE.
32 !
34 module fluid_pnpn
35  use num_types, only : rp
36  use krylov, only : ksp_monitor_t
38  pnpn_prs_res_factory, pnpn_vel_res_factory, &
39  pnpn_prs_res_stress_factory, pnpn_vel_res_stress_factory
41  rhs_maker_oifs_t, rhs_maker_sumab_fctry, rhs_maker_bdf_fctry, &
42  rhs_maker_ext_fctry, rhs_maker_oifs_fctry
43  use fluid_volflow, only : fluid_volflow_t
44  use fluid_scheme, only : fluid_scheme_t
46  use fluid_aux, only : fluid_step_info
48  use projection, only : projection_t
50  use advection, only : advection_t, advection_factory
53  use json_module, only : json_file
54  use ax_product, only : ax_t, ax_helm_factory
55  use field, only : field_t
56  use dirichlet, only : dirichlet_t
57  use facet_normal, only : facet_normal_t
58  use non_normal, only : non_normal_t
59  use mesh, only : mesh_t
60  use user_intf, only : user_t
62  use gs_ops, only : gs_op_add
63  use neko_config, only : neko_bcknd_device
64  use math, only : col2, glsum
65  use mathops, only : opadd2cm, opcolv
68  use utils, only : neko_error
69  use field_math, only : field_add2, field_copy
70  implicit none
71  private
72 
73 
74  type, public, extends(fluid_scheme_t) :: fluid_pnpn_t
75  type(field_t) :: p_res, u_res, v_res, w_res
76 
77  type(field_t) :: dp, du, dv, dw
78 
79  ! Coupled Helmholz operator for velocity
80  class(ax_t), allocatable :: ax_vel
81  ! Helmholz operator for pressure
82  class(ax_t), allocatable :: ax_prs
83 
84  type(projection_t) :: proj_prs
85  type(projection_t) :: proj_u
86  type(projection_t) :: proj_v
87  type(projection_t) :: proj_w
88 
89  type(facet_normal_t) :: bc_prs_surface
90  type(facet_normal_t) :: bc_sym_surface
91  type(dirichlet_t) :: bc_vel_res
92  type(dirichlet_t) :: bc_field_dirichlet_p
93  type(dirichlet_t) :: bc_field_dirichlet_u
94  type(dirichlet_t) :: bc_field_dirichlet_v
95  type(dirichlet_t) :: bc_field_dirichlet_w
96  type(non_normal_t) :: bc_vel_res_non_normal
97  type(bc_list_t) :: bclst_vel_res
98  type(bc_list_t) :: bclst_du
99  type(bc_list_t) :: bclst_dv
100  type(bc_list_t) :: bclst_dw
101  type(bc_list_t) :: bclst_dp
102 
103  class(advection_t), allocatable :: adv
104 
105  ! Time interpolation scheme
106  logical :: oifs
107 
108  ! Time variables
109  type(field_t) :: abx1, aby1, abz1
110  type(field_t) :: abx2, aby2, abz2
111  ! Advection terms for the oifs method
112  type(field_t) :: advx, advy, advz
113 
115  class(pnpn_prs_res_t), allocatable :: prs_res
116 
118  class(pnpn_vel_res_t), allocatable :: vel_res
119 
121  class(rhs_maker_sumab_t), allocatable :: sumab
122 
124  class(rhs_maker_ext_t), allocatable :: makeabf
125 
127  class(rhs_maker_bdf_t), allocatable :: makebdf
128 
130  class(rhs_maker_oifs_t), allocatable :: makeoifs
131 
133  type(fluid_volflow_t) :: vol_flow
134 
135  contains
136  procedure, pass(this) :: init => fluid_pnpn_init
137  procedure, pass(this) :: free => fluid_pnpn_free
138  procedure, pass(this) :: step => fluid_pnpn_step
139  procedure, pass(this) :: restart => fluid_pnpn_restart
140  end type fluid_pnpn_t
141 
142 contains
143 
144  subroutine fluid_pnpn_init(this, msh, lx, params, user, time_scheme)
145  class(fluid_pnpn_t), target, intent(inout) :: this
146  type(mesh_t), target, intent(inout) :: msh
147  integer, intent(inout) :: lx
148  type(json_file), target, intent(inout) :: params
149  type(user_t), target, intent(in) :: user
150  type(time_scheme_controller_t), target, intent(in) :: time_scheme
151  character(len=15), parameter :: scheme = 'Modular (Pn/Pn)'
152  logical :: advection
153 
154  call this%free()
155 
156  ! Initialize base class
157  call this%scheme_init(msh, lx, params, .true., .true., scheme, user)
158 
159  if (this%variable_material_properties .eqv. .true.) then
160  ! Setup backend dependent Ax routines
161  call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
162 
163  ! Setup backend dependent prs residual routines
164  call pnpn_prs_res_stress_factory(this%prs_res)
165 
166  ! Setup backend dependent vel residual routines
167  call pnpn_vel_res_stress_factory(this%vel_res)
168  else
169  ! Setup backend dependent Ax routines
170  call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
171 
172  ! Setup backend dependent prs residual routines
173  call pnpn_prs_res_factory(this%prs_res)
174 
175  ! Setup backend dependent vel residual routines
176  call pnpn_vel_res_factory(this%vel_res)
177  end if
178 
179  ! Setup Ax for the pressure
180  call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
181 
182 
183  ! Setup backend dependent summation of AB/BDF
184  call rhs_maker_sumab_fctry(this%sumab)
185 
186  ! Setup backend dependent summation of extrapolation scheme
187  call rhs_maker_ext_fctry(this%makeabf)
188 
189  ! Setup backend depenent contributions to F from lagged BD terms
190  call rhs_maker_bdf_fctry(this%makebdf)
191 
192  ! Setup backend dependent summations of the OIFS method
193  call rhs_maker_oifs_fctry(this%makeoifs)
194 
195  ! Initialize variables specific to this plan
196  associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
197  dm_xh => this%dm_Xh, nelv => this%msh%nelv)
198 
199  call this%p_res%init(dm_xh, "p_res")
200  call this%u_res%init(dm_xh, "u_res")
201  call this%v_res%init(dm_xh, "v_res")
202  call this%w_res%init(dm_xh, "w_res")
203  call this%abx1%init(dm_xh, "abx1")
204  call this%aby1%init(dm_xh, "aby1")
205  call this%abz1%init(dm_xh, "abz1")
206  call this%abx2%init(dm_xh, "abx2")
207  call this%aby2%init(dm_xh, "aby2")
208  call this%abz2%init(dm_xh, "abz2")
209  call this%advx%init(dm_xh, "advx")
210  call this%advy%init(dm_xh, "advy")
211  call this%advz%init(dm_xh, "advz")
212  this%abx1 = 0.0_rp
213  this%aby1 = 0.0_rp
214  this%abz1 = 0.0_rp
215  this%abx2 = 0.0_rp
216  this%aby2 = 0.0_rp
217  this%abz2 = 0.0_rp
218  this%advx = 0.0_rp
219  this%advy = 0.0_rp
220  this%advz = 0.0_rp
221 
222  call this%du%init(dm_xh, 'du')
223  call this%dv%init(dm_xh, 'dv')
224  call this%dw%init(dm_xh, 'dw')
225  call this%dp%init(dm_xh, 'dp')
226 
227  end associate
228 
229  ! Initialize velocity surface terms in pressure rhs
230  call this%bc_prs_surface%init_base(this%c_Xh)
231  call this%bc_prs_surface%mark_zone(msh%inlet)
232  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
233  'v', this%bc_labels)
234  ! This impacts the rhs of the pressure,
235  ! need to check what is correct to add here
236  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
237  'd_vel_u', this%bc_labels)
238  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
239  'd_vel_v', this%bc_labels)
240  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
241  'd_vel_w', this%bc_labels)
242  call this%bc_prs_surface%finalize()
243  ! Initialize symmetry surface terms in pressure rhs
244  call this%bc_sym_surface%init_base(this%c_Xh)
245  call this%bc_sym_surface%mark_zone(msh%sympln)
246  call this%bc_sym_surface%mark_zones_from_list(msh%labeled_zones,&
247  'sym', this%bc_labels)
248  ! Same here, should du, dv, dw be marked here?
249  call this%bc_sym_surface%finalize()
250  ! Initialize dirichlet bcs for velocity residual
251  call this%bc_vel_res_non_normal%init_base(this%c_Xh)
252  call this%bc_vel_res_non_normal%mark_zone(msh%outlet_normal)
253  call this%bc_vel_res_non_normal%mark_zones_from_list(msh%labeled_zones,&
254  'on', this%bc_labels)
255  call this%bc_vel_res_non_normal%mark_zones_from_list(msh%labeled_zones,&
256  'on+dong', &
257  this%bc_labels)
258  call this%bc_vel_res_non_normal%finalize()
259  call this%bc_vel_res_non_normal%init(this%c_Xh)
260 
261  call this%bc_field_dirichlet_p%init_base(this%c_Xh)
262  call this%bc_field_dirichlet_p%mark_zones_from_list(msh%labeled_zones, &
263  'on+dong', this%bc_labels)
264  call this%bc_field_dirichlet_p%mark_zones_from_list(msh%labeled_zones, &
265  'o+dong', this%bc_labels)
266  call this%bc_field_dirichlet_p%mark_zones_from_list(msh%labeled_zones, &
267  'd_pres', this%bc_labels)
268  call this%bc_field_dirichlet_p%finalize()
269  call this%bc_field_dirichlet_p%set_g(0.0_rp)
270  call bc_list_init(this%bclst_dp)
271  call bc_list_add(this%bclst_dp, this%bc_field_dirichlet_p)
272  !Add 0 prs bcs
273  call bc_list_add(this%bclst_dp, this%bc_prs)
274 
275  call this%bc_field_dirichlet_u%init_base(this%c_Xh)
276  call this%bc_field_dirichlet_u%mark_zones_from_list( &
277  msh%labeled_zones, 'd_vel_u', this%bc_labels)
278  call this%bc_field_dirichlet_u%finalize()
279  call this%bc_field_dirichlet_u%set_g(0.0_rp)
280 
281  call this%bc_field_dirichlet_v%init_base(this%c_Xh)
282  call this%bc_field_dirichlet_v%mark_zones_from_list(msh%labeled_zones, &
283  'd_vel_v', &
284  this%bc_labels)
285  call this%bc_field_dirichlet_v%finalize()
286  call this%bc_field_dirichlet_v%set_g(0.0_rp)
287 
288  call this%bc_field_dirichlet_w%init_base(this%c_Xh)
289  call this%bc_field_dirichlet_w%mark_zones_from_list(msh%labeled_zones, &
290  'd_vel_w', &
291  this%bc_labels)
292  call this%bc_field_dirichlet_w%finalize()
293  call this%bc_field_dirichlet_w%set_g(0.0_rp)
294 
295  call this%bc_vel_res%init_base(this%c_Xh)
296  call this%bc_vel_res%mark_zone(msh%inlet)
297  call this%bc_vel_res%mark_zone(msh%wall)
298  call this%bc_vel_res%mark_zones_from_list(msh%labeled_zones, &
299  'v', this%bc_labels)
300  call this%bc_vel_res%mark_zones_from_list(msh%labeled_zones, &
301  'w', this%bc_labels)
302  call this%bc_vel_res%finalize()
303  call this%bc_vel_res%set_g(0.0_rp)
304  call bc_list_init(this%bclst_vel_res)
305  call bc_list_add(this%bclst_vel_res, this%bc_vel_res)
306  call bc_list_add(this%bclst_vel_res, this%bc_vel_res_non_normal)
307  call bc_list_add(this%bclst_vel_res, this%bc_sym)
308  call bc_list_add(this%bclst_vel_res, this%bc_sh%symmetry)
309  call bc_list_add(this%bclst_vel_res, this%bc_wallmodel%symmetry)
310 
311  !Initialize bcs for u, v, w velocity components
312  call bc_list_init(this%bclst_du)
313  call bc_list_add(this%bclst_du, this%bc_sym%bc_x)
314  call bc_list_add(this%bclst_du, this%bc_sh%symmetry%bc_x)
315  call bc_list_add(this%bclst_du, this%bc_wallmodel%symmetry%bc_x)
316  call bc_list_add(this%bclst_du, this%bc_vel_res_non_normal%bc_x)
317  call bc_list_add(this%bclst_du, this%bc_vel_res)
318  call bc_list_add(this%bclst_du, this%bc_field_dirichlet_u)
319 
320  call bc_list_init(this%bclst_dv)
321  call bc_list_add(this%bclst_dv, this%bc_sym%bc_y)
322  call bc_list_add(this%bclst_dv, this%bc_sh%symmetry%bc_y)
323  call bc_list_add(this%bclst_dv, this%bc_wallmodel%symmetry%bc_y)
324  call bc_list_add(this%bclst_dv, this%bc_vel_res_non_normal%bc_y)
325  call bc_list_add(this%bclst_dv, this%bc_vel_res)
326  call bc_list_add(this%bclst_dv, this%bc_field_dirichlet_v)
327 
328  call bc_list_init(this%bclst_dw)
329  call bc_list_add(this%bclst_dw, this%bc_sym%bc_z)
330  call bc_list_add(this%bclst_dw, this%bc_sh%symmetry%bc_z)
331  call bc_list_add(this%bclst_dw, this%bc_wallmodel%symmetry%bc_z)
332  call bc_list_add(this%bclst_dw, this%bc_vel_res_non_normal%bc_z)
333  call bc_list_add(this%bclst_dw, this%bc_vel_res)
334  call bc_list_add(this%bclst_dw, this%bc_field_dirichlet_w)
335 
336  !Intialize projection space thingy
337 
338  if (this%variable_material_properties .and. &
339  this%vel_projection_dim .gt. 0) then
340  call neko_error("Velocity projection not available for full stress &
341  &formulation")
342  end if
343 
344 
345  call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
346  this%pr_projection_activ_step)
347 
348  call this%proj_u%init(this%dm_Xh%size(), this%vel_projection_dim, &
349  this%vel_projection_activ_step)
350  call this%proj_v%init(this%dm_Xh%size(), this%vel_projection_dim, &
351  this%vel_projection_activ_step)
352  call this%proj_w%init(this%dm_Xh%size(), this%vel_projection_dim, &
353  this%vel_projection_activ_step)
354 
355 
356  ! Add lagged term to checkpoint
357  call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
358 
359  ! Determine the time-interpolation scheme
360  call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
361 
362  ! Initialize the advection factory
363  call json_get_or_default(params, 'case.fluid.advection', advection, .true.)
364  call advection_factory(this%adv, params, this%c_Xh, &
365  this%ulag, this%vlag, this%wlag, &
366  this%chkp%dtlag, this%chkp%tlag, time_scheme, &
367  .not. advection)
368 
369  if (params%valid_path('case.fluid.flow_rate_force')) then
370  call this%vol_flow%init(this%dm_Xh, params)
371  end if
372 
373  end subroutine fluid_pnpn_init
374 
375  subroutine fluid_pnpn_restart(this, dtlag, tlag)
376  class(fluid_pnpn_t), target, intent(inout) :: this
377  real(kind=rp) :: dtlag(10), tlag(10)
378  type(field_t) :: u_temp, v_temp, w_temp
379  integer :: i, n
380 
381  n = this%u%dof%size()
382  if (allocated(this%chkp%previous_mesh%elements) .or. &
383  this%chkp%previous_Xh%lx .ne. this%Xh%lx) then
384  call col2(this%u%x, this%c_Xh%mult, this%u%dof%size())
385  call col2(this%v%x, this%c_Xh%mult, this%u%dof%size())
386  call col2(this%w%x, this%c_Xh%mult, this%u%dof%size())
387  call col2(this%p%x, this%c_Xh%mult, this%u%dof%size())
388  do i = 1, this%ulag%size()
389  call col2(this%ulag%lf(i)%x, this%c_Xh%mult, this%u%dof%size())
390  call col2(this%vlag%lf(i)%x, this%c_Xh%mult, this%u%dof%size())
391  call col2(this%wlag%lf(i)%x, this%c_Xh%mult, this%u%dof%size())
392  end do
393  end if
394 
395  if (neko_bcknd_device .eq. 1) then
396  associate(u => this%u, v => this%v, w => this%w, &
397  ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
398  p => this%p)
399  call device_memcpy(u%x, u%x_d, u%dof%size(), &
400  host_to_device, sync = .false.)
401  call device_memcpy(v%x, v%x_d, v%dof%size(), &
402  host_to_device, sync = .false.)
403  call device_memcpy(w%x, w%x_d, w%dof%size(), &
404  host_to_device, sync = .false.)
405  call device_memcpy(p%x, p%x_d, p%dof%size(), &
406  host_to_device, sync = .false.)
407  call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
408  u%dof%size(), host_to_device, sync = .false.)
409  call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
410  u%dof%size(), host_to_device, sync = .false.)
411 
412  call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
413  v%dof%size(), host_to_device, sync = .false.)
414  call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
415  v%dof%size(), host_to_device, sync = .false.)
416 
417  call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
418  w%dof%size(), host_to_device, sync = .false.)
419  call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
420  w%dof%size(), host_to_device, sync = .false.)
421  call device_memcpy(this%abx1%x, this%abx1%x_d, &
422  w%dof%size(), host_to_device, sync = .false.)
423  call device_memcpy(this%abx2%x, this%abx2%x_d, &
424  w%dof%size(), host_to_device, sync = .false.)
425  call device_memcpy(this%aby1%x, this%aby1%x_d, &
426  w%dof%size(), host_to_device, sync = .false.)
427  call device_memcpy(this%aby2%x, this%aby2%x_d, &
428  w%dof%size(), host_to_device, sync = .false.)
429  call device_memcpy(this%abz1%x, this%abz1%x_d, &
430  w%dof%size(), host_to_device, sync = .false.)
431  call device_memcpy(this%abz2%x, this%abz2%x_d, &
432  w%dof%size(), host_to_device, sync = .false.)
433  call device_memcpy(this%advx%x, this%advx%x_d, &
434  w%dof%size(), host_to_device, sync = .false.)
435  call device_memcpy(this%advy%x, this%advy%x_d, &
436  w%dof%size(), host_to_device, sync = .false.)
437  call device_memcpy(this%advz%x, this%advz%x_d, &
438  w%dof%size(), host_to_device, sync = .false.)
439  end associate
440  end if
441  ! Make sure that continuity is maintained (important for interpolation)
442  ! Do not do this for lagged rhs
443  ! (derivatives are not necessairly coninous across elements)
444 
445  if (allocated(this%chkp%previous_mesh%elements) &
446  .or. this%chkp%previous_Xh%lx .ne. this%Xh%lx) then
447  call this%gs_Xh%op(this%u, gs_op_add)
448  call this%gs_Xh%op(this%v, gs_op_add)
449  call this%gs_Xh%op(this%w, gs_op_add)
450  call this%gs_Xh%op(this%p, gs_op_add)
451 
452  do i = 1, this%ulag%size()
453  call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
454  call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
455  call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
456  end do
457  end if
458  !! If we would decide to only restart from lagged fields instead of saving
459  !! abx1, aby1 etc.
460  !! Observe that one also needs to recompute the focing at the old time steps
461  !u_temp = this%ulag%lf(2)
462  !v_temp = this%vlag%lf(2)
463  !w_temp = this%wlag%lf(2)
464  !! Compute the source terms
465  !call this%source_term%compute(tlag(2), -1)
466  !
467  !! Pre-multiply the source terms with the mass matrix.
468  !if (NEKO_BCKND_DEVICE .eq. 1) then
469  ! call device_opcolv(this%f_x%x_d, this%f_y%x_d, this%f_z%x_d, &
470  ! this%c_Xh%B_d, this%msh%gdim, n)
471  !else
472  ! call opcolv(this%f_x%x, this%f_y%x, this%f_z%x, &
473  ! this%c_Xh%B, this%msh%gdim, n)
474  !end if
475 
476  !! Add the advection operators to the right-hand-side.
477  !call this%adv%compute(u_temp, v_temp, w_temp, &
478  ! this%f_x%x, this%f_y%x, this%f_z%x, &
479  ! this%Xh, this%c_Xh, this%dm_Xh%size())
480  !this%abx2 = this%f_x
481  !this%aby2 = this%f_y
482  !this%abz2 = this%f_z
483  !
484  !u_temp = this%ulag%lf(1)
485  !v_temp = this%vlag%lf(1)
486  !w_temp = this%wlag%lf(1)
487  !call this%source_term%compute(tlag(1), 0)
488 
489  !! Pre-multiply the source terms with the mass matrix.
490  !if (NEKO_BCKND_DEVICE .eq. 1) then
491  ! call device_opcolv(this%f_x%x_d, this%f_y%x_d, this%f_z%x_d, &
492  ! this%c_Xh%B_d, this%msh%gdim, n)
493  !else
494  ! call opcolv(this%f_x%x, this%f_y%x, this%f_z%x, &
495  ! this%c_Xh%B, this%msh%gdim, n)
496  !end if
497 
498  !! Pre-multiply the source terms with the mass matrix.
499  !if (NEKO_BCKND_DEVICE .eq. 1) then
500  ! call device_opcolv(this%f_x%x_d, this%f_y%x_d, this%f_z%x_d, &
501  ! this%c_Xh%B_d, this%msh%gdim, n)
502  !else
503  ! call opcolv(this%f_x%x, this%f_y%x, this%f_z%x, &
504  ! this%c_Xh%B, this%msh%gdim, n)
505  !end if
506 
507  !call this%adv%compute(u_temp, v_temp, w_temp, &
508  ! this%f_x%x, this%f_y%x, this%f_z%x, &
509  ! this%Xh, this%c_Xh, this%dm_Xh%size())
510  !this%abx1 = this%f_x
511  !this%aby1 = this%f_y
512  !this%abz1 = this%f_z
513 
514  end subroutine fluid_pnpn_restart
515 
516  subroutine fluid_pnpn_free(this)
517  class(fluid_pnpn_t), intent(inout) :: this
518 
519  !Deallocate velocity and pressure fields
520  call this%scheme_free()
521 
522  call this%bc_prs_surface%free()
523  call this%bc_sym_surface%free()
524  call bc_list_free(this%bclst_vel_res)
525  call bc_list_free(this%bclst_dp)
526  call this%proj_prs%free()
527  call this%proj_u%free()
528  call this%proj_v%free()
529  call this%proj_w%free()
530 
531  call this%p_res%free()
532  call this%u_res%free()
533  call this%v_res%free()
534  call this%w_res%free()
535 
536  call this%du%free()
537  call this%dv%free()
538  call this%dw%free()
539  call this%dp%free()
540 
541  call this%abx1%free()
542  call this%aby1%free()
543  call this%abz1%free()
544 
545  call this%abx2%free()
546  call this%aby2%free()
547  call this%abz2%free()
548 
549  call this%advx%free()
550  call this%advy%free()
551  call this%advz%free()
552 
553  if (allocated(this%Ax_vel)) then
554  deallocate(this%Ax_vel)
555  end if
556 
557  if (allocated(this%Ax_prs)) then
558  deallocate(this%Ax_prs)
559  end if
560 
561  if (allocated(this%prs_res)) then
562  deallocate(this%prs_res)
563  end if
564 
565  if (allocated(this%vel_res)) then
566  deallocate(this%vel_res)
567  end if
568 
569  if (allocated(this%sumab)) then
570  deallocate(this%sumab)
571  end if
572 
573  if (allocated(this%makeabf)) then
574  deallocate(this%makeabf)
575  end if
576 
577  if (allocated(this%makebdf)) then
578  deallocate(this%makebdf)
579  end if
580 
581  if (allocated(this%makeoifs)) then
582  deallocate(this%makeoifs)
583  end if
584 
585  call this%vol_flow%free()
586 
587  end subroutine fluid_pnpn_free
588 
595  subroutine fluid_pnpn_step(this, t, tstep, dt, ext_bdf, dt_controller)
596  class(fluid_pnpn_t), target, intent(inout) :: this
597  real(kind=rp), intent(inout) :: t
598  integer, intent(inout) :: tstep
599  real(kind=rp), intent(in) :: dt
600  type(time_scheme_controller_t), intent(inout) :: ext_bdf
601  type(time_step_controller_t), intent(in) :: dt_controller
602  ! number of degrees of freedom
603  integer :: n
604  ! Solver results monitors (pressure + 3 velocity)
605  type(ksp_monitor_t) :: ksp_results(4)
606  ! Extrapolated velocity for the pressure residual
607  type(field_t), pointer :: u_e, v_e, w_e
608  ! Indices for tracking temporary fields
609  integer :: temp_indices(3)
610 
611  if (this%freeze) return
612 
613  n = this%dm_Xh%size()
614 
615  call profiler_start_region('Fluid', 1)
616  associate(u => this%u, v => this%v, w => this%w, p => this%p, &
617  du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
618  u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
619  p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
620  xh => this%Xh, &
621  c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
622  ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
623  msh => this%msh, prs_res => this%prs_res, &
624  source_term => this%source_term, vel_res => this%vel_res, &
625  sumab => this%sumab, makeoifs => this%makeoifs, &
626  makeabf => this%makeabf, makebdf => this%makebdf, &
627  vel_projection_dim => this%vel_projection_dim, &
628  pr_projection_dim => this%pr_projection_dim, &
629  rho => this%rho, mu => this%mu, oifs => this%oifs, &
630  rho_field => this%rho_field, mu_field => this%mu_field, &
631  f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
632  if_variable_dt => dt_controller%if_variable_dt, &
633  dt_last_change => dt_controller%dt_last_change)
634 
635  ! Get temporary arrays
636  call this%scratch%request_field(u_e, temp_indices(1))
637  call this%scratch%request_field(v_e, temp_indices(2))
638  call this%scratch%request_field(w_e, temp_indices(3))
639  call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
640  ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
641 
642  ! Compute the source terms
643  call this%source_term%compute(t, tstep)
644 
645  ! Add Neumann bc contributions to the RHS
646  call bc_list_apply_vector(this%bclst_vel_neumann, f_x%x, f_y%x, f_z%x, &
647  this%dm_Xh%size(), t, tstep)
648 
649  ! Compute the grandient jump penalty term
650  if (this%if_gradient_jump_penalty .eqv. .true.) then
651  call this%gradient_jump_penalty_u%compute(u, v, w, u)
652  call this%gradient_jump_penalty_v%compute(u, v, w, v)
653  call this%gradient_jump_penalty_w%compute(u, v, w, w)
654  call this%gradient_jump_penalty_u%perform(f_x)
655  call this%gradient_jump_penalty_v%perform(f_y)
656  call this%gradient_jump_penalty_w%perform(f_z)
657  end if
658 
659  if (oifs) then
660  ! Add the advection operators to the right-hand-side.
661  call this%adv%compute(u, v, w, &
662  this%advx, this%advy, this%advz, &
663  xh, this%c_Xh, dm_xh%size(), dt)
664 
665  ! At this point the RHS contains the sum of the advection operator and
666  ! additional source terms, evaluated using the velocity field from the
667  ! previous time-step. Now, this value is used in the explicit time
668  ! scheme to advance both terms in time.
669  call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
670  this%abx2, this%aby2, this%abz2, &
671  f_x%x, f_y%x, f_z%x, &
672  rho, ext_bdf%advection_coeffs, n)
673 
674  ! Now, the source terms from the previous time step are added to the RHS.
675  call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
676  f_x%x, f_y%x, f_z%x, &
677  rho, dt, n)
678  else
679  ! Add the advection operators to the right-hand-side.
680  call this%adv%compute(u, v, w, &
681  f_x, f_y, f_z, &
682  xh, this%c_Xh, dm_xh%size())
683 
684  ! At this point the RHS contains the sum of the advection operator and
685  ! additional source terms, evaluated using the velocity field from the
686  ! previous time-step. Now, this value is used in the explicit time
687  ! scheme to advance both terms in time.
688  call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
689  this%abx2, this%aby2, this%abz2, &
690  f_x%x, f_y%x, f_z%x, &
691  rho, ext_bdf%advection_coeffs, n)
692 
693  ! Add the RHS contributions coming from the BDF scheme.
694  call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
695  u, v, w, c_xh%B, rho, dt, &
696  ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
697  end if
698 
699  call ulag%update()
700  call vlag%update()
701  call wlag%update()
702 
706  call this%user_field_bc_vel%update(this%user_field_bc_vel%field_list, &
707  this%user_field_bc_vel%bc_list, this%c_Xh, t, tstep, "fluid")
708 
709  call this%bc_apply_vel(t, tstep)
710  call this%bc_apply_prs(t, tstep)
711 
712  ! Update material properties if necessary
713  call this%update_material_properties()
714 
715  ! Compute pressure.
716  call profiler_start_region('Pressure_residual', 18)
717  call prs_res%compute(p, p_res,&
718  u, v, w, &
719  u_e, v_e, w_e, &
720  f_x, f_y, f_z, &
721  c_xh, gs_xh, &
722  this%bc_prs_surface, this%bc_sym_surface,&
723  ax_prs, ext_bdf%diffusion_coeffs(1), dt, &
724  mu_field, rho_field)
725 
726  call gs_xh%op(p_res, gs_op_add)
727  call bc_list_apply_scalar(this%bclst_dp, p_res%x, p%dof%size(), t, tstep)
728  call profiler_end_region('Pressure_residual', 18)
729 
730  call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
731  'Pressure')
732 
733  call this%pc_prs%update()
734  call profiler_start_region('Pressure_solve', 3)
735  ksp_results(1) = &
736  this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, this%bclst_dp, gs_xh)
737 
738  call profiler_end_region('Pressure_solve', 3)
739 
740  call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
741  this%bclst_dp, gs_xh, n, tstep, dt_controller)
742 
743  call field_add2(p, dp, n)
744 
745  ! Compute velocity.
746  call profiler_start_region('Velocity_residual', 19)
747  call vel_res%compute(ax_vel, u, v, w, &
748  u_res, v_res, w_res, &
749  p, &
750  f_x, f_y, f_z, &
751  c_xh, msh, xh, &
752  mu_field, rho_field, ext_bdf%diffusion_coeffs(1), &
753  dt, dm_xh%size())
754 
755  call gs_xh%op(u_res, gs_op_add)
756  call gs_xh%op(v_res, gs_op_add)
757  call gs_xh%op(w_res, gs_op_add)
758 
759  call bc_list_apply_vector(this%bclst_vel_res,&
760  u_res%x, v_res%x, w_res%x, dm_xh%size(),&
761  t, tstep)
762 
763  ! We should implement a bc that takes three field_bcs and implements
764  ! vector_apply
765  if (neko_bcknd_device .eq. 1) then
766  call this%bc_field_dirichlet_u%apply_scalar_dev(u_res%x_d, t, tstep)
767  call this%bc_field_dirichlet_v%apply_scalar_dev(v_res%x_d, t, tstep)
768  call this%bc_field_dirichlet_w%apply_scalar_dev(w_res%x_d, t, tstep)
769  else
770  call this%bc_field_dirichlet_u%apply_scalar(u_res%x, n, t, tstep)
771  call this%bc_field_dirichlet_v%apply_scalar(v_res%x, n, t, tstep)
772  call this%bc_field_dirichlet_w%apply_scalar(w_res%x, n, t, tstep)
773  end if
774 
775  call profiler_end_region('Velocity_residual', 19)
776 
777  call this%proj_u%pre_solving(u_res%x, tstep, c_xh, n, dt_controller)
778  call this%proj_v%pre_solving(v_res%x, tstep, c_xh, n, dt_controller)
779  call this%proj_w%pre_solving(w_res%x, tstep, c_xh, n, dt_controller)
780 
781  call this%pc_vel%update()
782 
783  call profiler_start_region("Velocity_solve", 4)
784  ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
785  u_res%x, v_res%x, w_res%x, n, c_xh, &
786  this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
787  this%ksp_vel%max_iter)
788  call profiler_end_region("Velocity_solve", 4)
789 
790  call this%proj_u%post_solving(du%x, ax_vel, c_xh, &
791  this%bclst_du, gs_xh, n, tstep, dt_controller)
792  call this%proj_v%post_solving(dv%x, ax_vel, c_xh, &
793  this%bclst_dv, gs_xh, n, tstep, dt_controller)
794  call this%proj_w%post_solving(dw%x, ax_vel, c_xh, &
795  this%bclst_dw, gs_xh, n, tstep, dt_controller)
796 
797  if (neko_bcknd_device .eq. 1) then
798  call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
799  du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
800  else
801  call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
802  end if
803 
804  if (this%forced_flow_rate) then
805  call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
806  c_xh, gs_xh, ext_bdf, rho, mu, dt, &
807  this%bclst_dp, this%bclst_du, this%bclst_dv, &
808  this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
809  this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
810  this%ksp_vel%max_iter)
811  end if
812 
813  call fluid_step_info(tstep, t, dt, ksp_results)
814 
815  call this%scratch%relinquish_field(temp_indices)
816 
817  end associate
818  call profiler_end_region('Fluid', 1)
819  end subroutine fluid_pnpn_step
820 
821 
822 end module fluid_pnpn
Copy data between host and device (or device and device)
Definition: device.F90:51
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Definition: json_utils.f90:54
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:45
Subroutines to add advection terms to the RHS of a transport equation.
Definition: advection.f90:34
Defines a Matrix-vector product.
Definition: ax.f90:34
Defines a boundary condition.
Definition: bc.f90:34
subroutine, public bc_list_add(bclst, bc)
Add a condition to a list of boundary conditions.
Definition: bc.f90:505
subroutine, public bc_list_apply_vector(bclst, x, y, z, n, t, tstep)
Apply a list of boundary conditions to a vector field.
Definition: bc.f90:587
subroutine, public bc_list_init(bclst, size)
Constructor for a list of boundary conditions.
Definition: bc.f90:464
subroutine, public bc_list_free(bclst)
Destructor for a list of boundary conditions.
Definition: bc.f90:491
subroutine, public bc_list_apply_scalar(bclst, x, n, t, tstep)
Apply a list of boundary conditions to a scalar field.
Definition: bc.f90:530
subroutine, public device_opadd2cm(a1_d, a2_d, a3_d, b1_d, b2_d, b3_d, c, n, gdim)
subroutine, public device_opcolv(a1_d, a2_d, a3_d, c_d, gdim, n)
Device abstraction, common interface for various accelerators.
Definition: device.F90:34
integer, parameter, public host_to_device
Definition: device.F90:47
Defines a dirichlet boundary condition.
Definition: dirichlet.f90:34
Dirichlet condition applied in the facet normal direction.
subroutine, public field_add2(a, b, n)
Vector addition .
Definition: field_math.f90:255
subroutine, public field_copy(a, b, n)
Copy a vector .
Definition: field_math.f90:126
Defines a field.
Definition: field.f90:34
Auxiliary routines for fluid solvers.
Definition: fluid_aux.f90:2
subroutine, public fluid_step_info(step, t, dt, ksp_results)
Prints for prs, velx, vely, velz the following: Number of iterations, start residual,...
Definition: fluid_aux.f90:17
Modular version of the Classic Nek5000 Pn/Pn formulation for fluids.
Definition: fluid_pnpn.f90:34
subroutine fluid_pnpn_init(this, msh, lx, params, user, time_scheme)
Definition: fluid_pnpn.f90:145
subroutine fluid_pnpn_step(this, t, tstep, dt, ext_bdf, dt_controller)
Advance fluid simulation in time.
Definition: fluid_pnpn.f90:596
subroutine fluid_pnpn_free(this)
Definition: fluid_pnpn.f90:517
subroutine fluid_pnpn_restart(this, dtlag, tlag)
Definition: fluid_pnpn.f90:376
Fluid formulations.
Defines Gather-scatter operations.
Definition: gs_ops.f90:34
integer, parameter, public gs_op_add
Definition: gs_ops.f90:36
Utilities for retrieving parameters from the case files.
Definition: json_utils.f90:34
Implements the base abstract type for Krylov solvers plus helper types.
Definition: krylov.f90:34
Definition: math.f90:60
real(kind=rp) function, public glsum(a, n)
Sum a vector of length n.
Definition: math.f90:360
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:729
Collection of vector field operations operating on and . Note that in general the indices and ....
Definition: mathops.f90:65
subroutine, public opcolv(a1, a2, a3, c, gdim, n)
Definition: mathops.f90:97
subroutine, public opadd2cm(a1, a2, a3, b1, b2, b3, c, n, gdim)
Definition: mathops.f90:142
Defines a mesh.
Definition: mesh.f90:34
Build configurations.
Definition: neko_config.f90:34
integer, parameter neko_bcknd_device
Definition: neko_config.f90:44
Dirichlet condition on axis aligned plane in the non normal direction.
Definition: non_normal.f90:34
integer, parameter, public rp
Global precision used in computations.
Definition: num_types.f90:12
Defines Pressure and velocity residuals in the Pn-Pn formulation.
Definition: pnpn_res.f90:34
Profiling interface.
Definition: profiler.F90:34
subroutine, public profiler_start_region(name, region_id)
Started a named (name) profiler region.
Definition: profiler.F90:78
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition: profiler.F90:109
Project x onto X, the space of old solutions and back again.
Definition: projection.f90:63
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition: rhs_maker.f90:38
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Definition: source_term.f90:34
Compound scheme for the advection and diffusion operators in a transport equation.
Base class for time integration schemes.
Definition: time_scheme.f90:61
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition: user_intf.f90:34
Utilities.
Definition: utils.f90:35
Base abstract type for computing the advection operator.
Definition: advection.f90:46
Base type for a matrix-vector product providing .
Definition: ax.f90:43
A list of boundary conditions.
Definition: bc.f90:104
Generic Dirichlet boundary condition on .
Definition: dirichlet.f90:44
Dirichlet condition in facet normal direction.
Base type of all fluid formulations.
Type for storing initial and final residuals in a Krylov solver.
Definition: krylov.f90:56
Dirichlet condition in non normal direction of a plane.
Definition: non_normal.f90:44
Abstract type to compute pressure residual.
Definition: pnpn_res.f90:47
Abstract type to compute velocity residual.
Definition: pnpn_res.f90:53
Abstract type to add contributions to F from lagged BD terms.
Definition: rhs_maker.f90:59
Abstract type to sum up contributions to kth order extrapolation scheme.
Definition: rhs_maker.f90:52
Abstract type to add contributions of kth order OIFS scheme.
Definition: rhs_maker.f90:66
Abstract type to compute extrapolated velocity field for the pressure equation.
Definition: rhs_maker.f90:46
Implements the logic to compute the time coefficients for the advection and diffusion operators in a ...
Provides a tool to set time step dt.