Neko  0.8.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
55  use ax_product, only : ax_t, ax_helm_factory
56  use field, only : field_t
57  use dirichlet, only : dirichlet_t
58  use facet_normal, only : facet_normal_t
59  use non_normal, only : non_normal_t
60  use mesh, only : mesh_t
61  use user_intf, only : user_t
63  use gs_ops, only : gs_op_add
64  use neko_config, only : neko_bcknd_device
65  use math, only : col2
66  use mathops, only : opadd2cm, opcolv
69  use utils, only : neko_error
70  use field_math, only : field_add2
71  implicit none
72  private
73 
74 
75  type, public, extends(fluid_scheme_t) :: fluid_pnpn_t
76  type(field_t) :: p_res, u_res, v_res, w_res
77 
78  type(field_t) :: dp, du, dv, dw
79 
80  ! Coupled Helmholz operator for velocity
81  class(ax_t), allocatable :: ax_vel
82  ! Helmholz operator for pressure
83  class(ax_t), allocatable :: ax_prs
84 
85  type(projection_t) :: proj_prs
86  type(projection_t) :: proj_u
87  type(projection_t) :: proj_v
88  type(projection_t) :: proj_w
89 
90  type(facet_normal_t) :: bc_prs_surface
91  type(facet_normal_t) :: bc_sym_surface
92  type(dirichlet_t) :: bc_vel_res
93  type(dirichlet_t) :: bc_field_dirichlet_p
94  type(dirichlet_t) :: bc_field_dirichlet_u
95  type(dirichlet_t) :: bc_field_dirichlet_v
96  type(dirichlet_t) :: bc_field_dirichlet_w
97  type(non_normal_t) :: bc_vel_res_non_normal
98  type(bc_list_t) :: bclst_vel_res
99  type(bc_list_t) :: bclst_du
100  type(bc_list_t) :: bclst_dv
101  type(bc_list_t) :: bclst_dw
102  type(bc_list_t) :: bclst_dp
103 
104  class(advection_t), allocatable :: adv
105 
106  ! Time interpolation scheme
107  logical :: oifs
108 
109  ! Time variables
110  type(field_t) :: abx1, aby1, abz1
111  type(field_t) :: abx2, aby2, abz2
112  ! Advection terms for the oifs method
113  type(field_t) :: advx, advy, advz
114 
116  class(pnpn_prs_res_t), allocatable :: prs_res
117 
119  class(pnpn_vel_res_t), allocatable :: vel_res
120 
122  class(rhs_maker_sumab_t), allocatable :: sumab
123 
125  class(rhs_maker_ext_t), allocatable :: makeabf
126 
128  class(rhs_maker_bdf_t), allocatable :: makebdf
129 
131  class(rhs_maker_oifs_t), allocatable :: makeoifs
132 
134  type(fluid_volflow_t) :: vol_flow
135 
136  contains
137  procedure, pass(this) :: init => fluid_pnpn_init
138  procedure, pass(this) :: free => fluid_pnpn_free
139  procedure, pass(this) :: step => fluid_pnpn_step
140  procedure, pass(this) :: restart => fluid_pnpn_restart
141  end type fluid_pnpn_t
142 
143 contains
144 
145  subroutine fluid_pnpn_init(this, msh, lx, params, user, &
146  material_properties, time_scheme)
147  class(fluid_pnpn_t), target, intent(inout) :: this
148  type(mesh_t), target, intent(inout) :: msh
149  integer, intent(inout) :: lx
150  type(json_file), target, intent(inout) :: params
151  type(user_t), intent(in) :: user
152  type(material_properties_t), target, intent(inout) :: material_properties
153  type(time_scheme_controller_t), target, intent(in) :: time_scheme
154  character(len=15), parameter :: scheme = 'Modular (Pn/Pn)'
155 
156  call this%free()
157 
158  ! Initialize base class
159  call this%scheme_init(msh, lx, params, .true., .true., scheme, user, &
161 
162  if (this%variable_material_properties .eqv. .true.) then
163  ! Setup backend dependent Ax routines
164  call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
165 
166  ! Setup backend dependent prs residual routines
167  call pnpn_prs_res_stress_factory(this%prs_res)
168 
169  ! Setup backend dependent vel residual routines
170  call pnpn_vel_res_stress_factory(this%vel_res)
171  else
172  ! Setup backend dependent Ax routines
173  call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
174 
175  ! Setup backend dependent prs residual routines
176  call pnpn_prs_res_factory(this%prs_res)
177 
178  ! Setup backend dependent vel residual routines
179  call pnpn_vel_res_factory(this%vel_res)
180  end if
181 
182  ! Setup Ax for the pressure
183  call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
184 
185 
186  ! Setup backend dependent summation of AB/BDF
187  call rhs_maker_sumab_fctry(this%sumab)
188 
189  ! Setup backend dependent summation of extrapolation scheme
190  call rhs_maker_ext_fctry(this%makeabf)
191 
192  ! Setup backend depenent contributions to F from lagged BD terms
193  call rhs_maker_bdf_fctry(this%makebdf)
194 
195  ! Setup backend dependent summations of the OIFS method
196  call rhs_maker_oifs_fctry(this%makeoifs)
197 
198  ! Initialize variables specific to this plan
199  associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
200  dm_xh => this%dm_Xh, nelv => this%msh%nelv)
201 
202  call this%p_res%init(dm_xh, "p_res")
203  call this%u_res%init(dm_xh, "u_res")
204  call this%v_res%init(dm_xh, "v_res")
205  call this%w_res%init(dm_xh, "w_res")
206  call this%abx1%init(dm_xh, "abx1")
207  call this%aby1%init(dm_xh, "aby1")
208  call this%abz1%init(dm_xh, "abz1")
209  call this%abx2%init(dm_xh, "abx2")
210  call this%aby2%init(dm_xh, "aby2")
211  call this%abz2%init(dm_xh, "abz2")
212  call this%advx%init(dm_xh, "advx")
213  call this%advy%init(dm_xh, "advy")
214  call this%advz%init(dm_xh, "advz")
215  this%abx1 = 0.0_rp
216  this%aby1 = 0.0_rp
217  this%abz1 = 0.0_rp
218  this%abx2 = 0.0_rp
219  this%aby2 = 0.0_rp
220  this%abz2 = 0.0_rp
221  this%advx = 0.0_rp
222  this%advy = 0.0_rp
223  this%advz = 0.0_rp
224 
225  call this%du%init(dm_xh, 'du')
226  call this%dv%init(dm_xh, 'dv')
227  call this%dw%init(dm_xh, 'dw')
228  call this%dp%init(dm_xh, 'dp')
229 
230  end associate
231 
232  ! Initialize velocity surface terms in pressure rhs
233  call this%bc_prs_surface%init_base(this%c_Xh)
234  call this%bc_prs_surface%mark_zone(msh%inlet)
235  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
236  'v', this%bc_labels)
237  ! This impacts the rhs of the pressure,
238  ! need to check what is correct to add here
239  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
240  'd_vel_u', this%bc_labels)
241  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
242  'd_vel_v', this%bc_labels)
243  call this%bc_prs_surface%mark_zones_from_list(msh%labeled_zones,&
244  'd_vel_w', this%bc_labels)
245  call this%bc_prs_surface%finalize()
246  ! Initialize symmetry surface terms in pressure rhs
247  call this%bc_sym_surface%init_base(this%c_Xh)
248  call this%bc_sym_surface%mark_zone(msh%sympln)
249  call this%bc_sym_surface%mark_zones_from_list(msh%labeled_zones,&
250  'sym', this%bc_labels)
251  ! Same here, should du, dv, dw be marked here?
252  call this%bc_sym_surface%finalize()
253  ! Initialize dirichlet bcs for velocity residual
254  call this%bc_vel_res_non_normal%init_base(this%c_Xh)
255  call this%bc_vel_res_non_normal%mark_zone(msh%outlet_normal)
256  call this%bc_vel_res_non_normal%mark_zones_from_list(msh%labeled_zones,&
257  'on', this%bc_labels)
258  call this%bc_vel_res_non_normal%mark_zones_from_list(msh%labeled_zones,&
259  'on+dong', &
260  this%bc_labels)
261  call this%bc_vel_res_non_normal%finalize()
262  call this%bc_vel_res_non_normal%init(this%c_Xh)
263 
264  call this%bc_field_dirichlet_p%init_base(this%c_Xh)
265  call this%bc_field_dirichlet_p%mark_zones_from_list(msh%labeled_zones, &
266  'on+dong', this%bc_labels)
267  call this%bc_field_dirichlet_p%mark_zones_from_list(msh%labeled_zones, &
268  'o+dong', this%bc_labels)
269  call this%bc_field_dirichlet_p%mark_zones_from_list(msh%labeled_zones, &
270  'd_pres', this%bc_labels)
271  call this%bc_field_dirichlet_p%finalize()
272  call this%bc_field_dirichlet_p%set_g(0.0_rp)
273  call bc_list_init(this%bclst_dp)
274  call bc_list_add(this%bclst_dp, this%bc_field_dirichlet_p)
275  !Add 0 prs bcs
276  call bc_list_add(this%bclst_dp, this%bc_prs)
277 
278  call this%bc_field_dirichlet_u%init_base(this%c_Xh)
279  call this%bc_field_dirichlet_u%mark_zones_from_list( &
280  msh%labeled_zones, 'd_vel_u', this%bc_labels)
281  call this%bc_field_dirichlet_u%finalize()
282  call this%bc_field_dirichlet_u%set_g(0.0_rp)
283 
284  call this%bc_field_dirichlet_v%init_base(this%c_Xh)
285  call this%bc_field_dirichlet_v%mark_zones_from_list(msh%labeled_zones, &
286  'd_vel_v', &
287  this%bc_labels)
288  call this%bc_field_dirichlet_v%finalize()
289  call this%bc_field_dirichlet_v%set_g(0.0_rp)
290 
291  call this%bc_field_dirichlet_w%init_base(this%c_Xh)
292  call this%bc_field_dirichlet_w%mark_zones_from_list(msh%labeled_zones, &
293  'd_vel_w', &
294  this%bc_labels)
295  call this%bc_field_dirichlet_w%finalize()
296  call this%bc_field_dirichlet_w%set_g(0.0_rp)
297 
298  call this%bc_vel_res%init_base(this%c_Xh)
299  call this%bc_vel_res%mark_zone(msh%inlet)
300  call this%bc_vel_res%mark_zone(msh%wall)
301  call this%bc_vel_res%mark_zones_from_list(msh%labeled_zones, &
302  'v', this%bc_labels)
303  call this%bc_vel_res%mark_zones_from_list(msh%labeled_zones, &
304  'w', this%bc_labels)
305  call this%bc_vel_res%finalize()
306  call this%bc_vel_res%set_g(0.0_rp)
307  call bc_list_init(this%bclst_vel_res)
308  call bc_list_add(this%bclst_vel_res, this%bc_vel_res)
309  call bc_list_add(this%bclst_vel_res, this%bc_vel_res_non_normal)
310  call bc_list_add(this%bclst_vel_res, this%bc_sym)
311 
312  !Initialize bcs for u, v, w velocity components
313  call bc_list_init(this%bclst_du)
314  call bc_list_add(this%bclst_du, this%bc_sym%bc_x)
315  call bc_list_add(this%bclst_du, this%bc_vel_res_non_normal%bc_x)
316  call bc_list_add(this%bclst_du, this%bc_vel_res)
317  call bc_list_add(this%bclst_du, this%bc_field_dirichlet_u)
318 
319  call bc_list_init(this%bclst_dv)
320  call bc_list_add(this%bclst_dv, this%bc_sym%bc_y)
321  call bc_list_add(this%bclst_dv, this%bc_vel_res_non_normal%bc_y)
322  call bc_list_add(this%bclst_dv, this%bc_vel_res)
323  call bc_list_add(this%bclst_dv, this%bc_field_dirichlet_v)
324 
325  call bc_list_init(this%bclst_dw)
326  call bc_list_add(this%bclst_dw, this%bc_sym%bc_z)
327  call bc_list_add(this%bclst_dw, this%bc_vel_res_non_normal%bc_z)
328  call bc_list_add(this%bclst_dw, this%bc_vel_res)
329  call bc_list_add(this%bclst_dw, this%bc_field_dirichlet_w)
330 
331  !Intialize projection space thingy
332 
333  if (this%variable_material_properties .and. &
334  this%vel_projection_dim .gt. 0) then
335  call neko_error("Velocity projection not available for full stress &
336  &formulation")
337  end if
338 
339 
340  call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
341  this%pr_projection_activ_step)
342 
343  call this%proj_u%init(this%dm_Xh%size(), this%vel_projection_dim, &
344  this%vel_projection_activ_step)
345  call this%proj_v%init(this%dm_Xh%size(), this%vel_projection_dim, &
346  this%vel_projection_activ_step)
347  call this%proj_w%init(this%dm_Xh%size(), this%vel_projection_dim, &
348  this%vel_projection_activ_step)
349 
350 
351  ! Add lagged term to checkpoint
352  call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
353 
354  ! Determine the time-interpolation scheme
355  call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
356 
357  ! Initialize the advection factory
358  call advection_factory(this%adv, params, this%c_Xh, &
359  this%ulag, this%vlag, this%wlag, &
360  this%chkp%dtlag, this%chkp%tlag, time_scheme)
361 
362  if (params%valid_path('case.fluid.flow_rate_force')) then
363  call this%vol_flow%init(this%dm_Xh, params)
364  end if
365 
366  end subroutine fluid_pnpn_init
367 
368  subroutine fluid_pnpn_restart(this, dtlag, tlag)
369  class(fluid_pnpn_t), target, intent(inout) :: this
370  real(kind=rp) :: dtlag(10), tlag(10)
371  type(field_t) :: u_temp, v_temp, w_temp
372  integer :: i, n
373 
374  n = this%u%dof%size()
375  ! Make sure that continuity is maintained (important for interpolation)
376  ! Do not do this for lagged rhs
377  ! (derivatives are not necessairly coninous across elements)
378  call col2(this%u%x, this%c_Xh%mult, this%u%dof%size())
379  call col2(this%v%x, this%c_Xh%mult, this%u%dof%size())
380  call col2(this%w%x, this%c_Xh%mult, this%u%dof%size())
381  call col2(this%p%x, this%c_Xh%mult, this%u%dof%size())
382  do i = 1, this%ulag%size()
383  call col2(this%ulag%lf(i)%x, this%c_Xh%mult, this%u%dof%size())
384  call col2(this%vlag%lf(i)%x, this%c_Xh%mult, this%u%dof%size())
385  call col2(this%wlag%lf(i)%x, this%c_Xh%mult, this%u%dof%size())
386  end do
387 
388  if (neko_bcknd_device .eq. 1) then
389  associate(u => this%u, v => this%v, w => this%w, &
390  ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
391  p => this%p)
392  call device_memcpy(u%x, u%x_d, u%dof%size(), &
393  host_to_device, sync = .false.)
394  call device_memcpy(v%x, v%x_d, v%dof%size(), &
395  host_to_device, sync = .false.)
396  call device_memcpy(w%x, w%x_d, w%dof%size(), &
397  host_to_device, sync = .false.)
398  call device_memcpy(p%x, p%x_d, p%dof%size(), &
399  host_to_device, sync = .false.)
400  call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
401  u%dof%size(), host_to_device, sync = .false.)
402  call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
403  u%dof%size(), host_to_device, sync = .false.)
404 
405  call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
406  v%dof%size(), host_to_device, sync = .false.)
407  call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
408  v%dof%size(), host_to_device, sync = .false.)
409 
410  call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
411  w%dof%size(), host_to_device, sync = .false.)
412  call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
413  w%dof%size(), host_to_device, sync = .false.)
414  call device_memcpy(this%abx1%x, this%abx1%x_d, &
415  w%dof%size(), host_to_device, sync = .false.)
416  call device_memcpy(this%abx2%x, this%abx2%x_d, &
417  w%dof%size(), host_to_device, sync = .false.)
418  call device_memcpy(this%aby1%x, this%aby1%x_d, &
419  w%dof%size(), host_to_device, sync = .false.)
420  call device_memcpy(this%aby2%x, this%aby2%x_d, &
421  w%dof%size(), host_to_device, sync = .false.)
422  call device_memcpy(this%abz1%x, this%abz1%x_d, &
423  w%dof%size(), host_to_device, sync = .false.)
424  call device_memcpy(this%abz2%x, this%abz2%x_d, &
425  w%dof%size(), host_to_device, sync = .false.)
426  call device_memcpy(this%advx%x, this%advx%x_d, &
427  w%dof%size(), host_to_device, sync = .false.)
428  call device_memcpy(this%advy%x, this%advy%x_d, &
429  w%dof%size(), host_to_device, sync = .false.)
430  call device_memcpy(this%advz%x, this%advz%x_d, &
431  w%dof%size(), host_to_device, sync = .false.)
432  end associate
433  end if
434 
435 
436  call this%gs_Xh%op(this%u, gs_op_add)
437  call this%gs_Xh%op(this%v, gs_op_add)
438  call this%gs_Xh%op(this%w, gs_op_add)
439  call this%gs_Xh%op(this%p, gs_op_add)
440 
441  do i = 1, this%ulag%size()
442  call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
443  call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
444  call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
445  end do
446 
447  !! If we would decide to only restart from lagged fields instead of saving
448  !! abx1, aby1 etc.
449  !! Observe that one also needs to recompute the focing at the old time steps
450  !u_temp = this%ulag%lf(2)
451  !v_temp = this%vlag%lf(2)
452  !w_temp = this%wlag%lf(2)
453  !! Compute the source terms
454  !call this%source_term%compute(tlag(2), -1)
455  !
456  !! Pre-multiply the source terms with the mass matrix.
457  !if (NEKO_BCKND_DEVICE .eq. 1) then
458  ! call device_opcolv(this%f_x%x_d, this%f_y%x_d, this%f_z%x_d, &
459  ! this%c_Xh%B_d, this%msh%gdim, n)
460  !else
461  ! call opcolv(this%f_x%x, this%f_y%x, this%f_z%x, &
462  ! this%c_Xh%B, this%msh%gdim, n)
463  !end if
464 
465  !! Add the advection operators to the right-hand-side.
466  !call this%adv%compute(u_temp, v_temp, w_temp, &
467  ! this%f_x%x, this%f_y%x, this%f_z%x, &
468  ! this%Xh, this%c_Xh, this%dm_Xh%size())
469  !this%abx2 = this%f_x
470  !this%aby2 = this%f_y
471  !this%abz2 = this%f_z
472  !
473  !u_temp = this%ulag%lf(1)
474  !v_temp = this%vlag%lf(1)
475  !w_temp = this%wlag%lf(1)
476  !call this%source_term%compute(tlag(1), 0)
477 
478  !! Pre-multiply the source terms with the mass matrix.
479  !if (NEKO_BCKND_DEVICE .eq. 1) then
480  ! call device_opcolv(this%f_x%x_d, this%f_y%x_d, this%f_z%x_d, &
481  ! this%c_Xh%B_d, this%msh%gdim, n)
482  !else
483  ! call opcolv(this%f_x%x, this%f_y%x, this%f_z%x, &
484  ! this%c_Xh%B, this%msh%gdim, n)
485  !end if
486 
487  !! Pre-multiply the source terms with the mass matrix.
488  !if (NEKO_BCKND_DEVICE .eq. 1) then
489  ! call device_opcolv(this%f_x%x_d, this%f_y%x_d, this%f_z%x_d, &
490  ! this%c_Xh%B_d, this%msh%gdim, n)
491  !else
492  ! call opcolv(this%f_x%x, this%f_y%x, this%f_z%x, &
493  ! this%c_Xh%B, this%msh%gdim, n)
494  !end if
495 
496  !call this%adv%compute(u_temp, v_temp, w_temp, &
497  ! this%f_x%x, this%f_y%x, this%f_z%x, &
498  ! this%Xh, this%c_Xh, this%dm_Xh%size())
499  !this%abx1 = this%f_x
500  !this%aby1 = this%f_y
501  !this%abz1 = this%f_z
502 
503  end subroutine fluid_pnpn_restart
504 
505  subroutine fluid_pnpn_free(this)
506  class(fluid_pnpn_t), intent(inout) :: this
507 
508  !Deallocate velocity and pressure fields
509  call this%scheme_free()
510 
511  call this%bc_prs_surface%free()
512  call this%bc_sym_surface%free()
513  call bc_list_free(this%bclst_vel_res)
514  call bc_list_free(this%bclst_dp)
515  call this%proj_prs%free()
516  call this%proj_u%free()
517  call this%proj_v%free()
518  call this%proj_w%free()
519 
520  call this%p_res%free()
521  call this%u_res%free()
522  call this%v_res%free()
523  call this%w_res%free()
524 
525  call this%du%free()
526  call this%dv%free()
527  call this%dw%free()
528  call this%dp%free()
529 
530  call this%abx1%free()
531  call this%aby1%free()
532  call this%abz1%free()
533 
534  call this%abx2%free()
535  call this%aby2%free()
536  call this%abz2%free()
537 
538  call this%advx%free()
539  call this%advy%free()
540  call this%advz%free()
541 
542  if (allocated(this%Ax_vel)) then
543  deallocate(this%Ax_vel)
544  end if
545 
546  if (allocated(this%Ax_prs)) then
547  deallocate(this%Ax_prs)
548  end if
549 
550  if (allocated(this%prs_res)) then
551  deallocate(this%prs_res)
552  end if
553 
554  if (allocated(this%vel_res)) then
555  deallocate(this%vel_res)
556  end if
557 
558  if (allocated(this%sumab)) then
559  deallocate(this%sumab)
560  end if
561 
562  if (allocated(this%makeabf)) then
563  deallocate(this%makeabf)
564  end if
565 
566  if (allocated(this%makebdf)) then
567  deallocate(this%makebdf)
568  end if
569 
570  if (allocated(this%makeoifs)) then
571  deallocate(this%makeoifs)
572  end if
573 
574  call this%vol_flow%free()
575 
576  end subroutine fluid_pnpn_free
577 
584  subroutine fluid_pnpn_step(this, t, tstep, dt, ext_bdf, dt_controller)
585  class(fluid_pnpn_t), target, intent(inout) :: this
586  real(kind=rp), intent(inout) :: t
587  integer, intent(inout) :: tstep
588  real(kind=rp), intent(in) :: dt
589  type(time_scheme_controller_t), intent(inout) :: ext_bdf
590  type(time_step_controller_t), intent(in) :: dt_controller
591  ! number of degrees of freedom
592  integer :: n
593  ! Solver results monitors (pressure + 3 velocity)
594  type(ksp_monitor_t) :: ksp_results(4)
595  ! Extrapolated velocity for the pressure residual
596  type(field_t), pointer :: u_e, v_e, w_e
597  ! Indices for tracking temporary fields
598  integer :: temp_indices(3)
599 
600  if (this%freeze) return
601 
602  n = this%dm_Xh%size()
603 
604  call profiler_start_region('Fluid', 1)
605  associate(u => this%u, v => this%v, w => this%w, p => this%p, &
606  du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
607  u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
608  p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
609  xh => this%Xh, &
610  c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
611  ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
612  msh => this%msh, prs_res => this%prs_res, &
613  source_term => this%source_term, vel_res => this%vel_res, &
614  sumab => this%sumab, makeoifs => this%makeoifs, &
615  makeabf => this%makeabf, makebdf => this%makebdf, &
616  vel_projection_dim => this%vel_projection_dim, &
617  pr_projection_dim => this%pr_projection_dim, &
618  rho => this%rho, mu => this%mu, oifs => this%oifs, &
619  rho_field => this%rho_field, mu_field => this%mu_field, &
620  f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
621  if_variable_dt => dt_controller%if_variable_dt, &
622  dt_last_change => dt_controller%dt_last_change)
623 
624  ! Get temporary arrays
625  call this%scratch%request_field(u_e, temp_indices(1))
626  call this%scratch%request_field(v_e, temp_indices(2))
627  call this%scratch%request_field(w_e, temp_indices(3))
628  call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
629  ulag, vlag, wlag, ext_bdf%advection_coeffs, ext_bdf%nadv)
630 
631  ! Compute the source terms
632  call this%source_term%compute(t, tstep)
633 
634  ! Pre-multiply the source terms with the mass matrix.
635  if (neko_bcknd_device .eq. 1) then
636  call device_opcolv(f_x%x_d, f_y%x_d, f_z%x_d, c_xh%B_d, msh%gdim, n)
637  else
638  call opcolv(f_x%x, f_y%x, f_z%x, c_xh%B, msh%gdim, n)
639  end if
640 
641  if (oifs) then
642  ! Add the advection operators to the right-hand-side.
643  call this%adv%compute(u, v, w, &
644  this%advx, this%advy, this%advz, &
645  xh, this%c_Xh, dm_xh%size(), dt)
646 
647  ! At this point the RHS contains the sum of the advection operator and
648  ! additional source terms, evaluated using the velocity field from the
649  ! previous time-step. Now, this value is used in the explicit time
650  ! scheme to advance both terms in time.
651  call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
652  this%abx2, this%aby2, this%abz2, &
653  f_x%x, f_y%x, f_z%x, &
654  rho, ext_bdf%advection_coeffs, n)
655 
656  ! Now, the source terms from the previous time step are added to the RHS.
657  call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
658  f_x%x, f_y%x, f_z%x, &
659  rho, dt, n)
660  else
661  ! Add the advection operators to the right-hand-side.
662  call this%adv%compute(u, v, w, &
663  f_x, f_y, f_z, &
664  xh, this%c_Xh, dm_xh%size())
665 
666  ! At this point the RHS contains the sum of the advection operator and
667  ! additional source terms, evaluated using the velocity field from the
668  ! previous time-step. Now, this value is used in the explicit time
669  ! scheme to advance both terms in time.
670  call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
671  this%abx2, this%aby2, this%abz2, &
672  f_x%x, f_y%x, f_z%x, &
673  rho, ext_bdf%advection_coeffs, n)
674 
675  ! Add the RHS contributions coming from the BDF scheme.
676  call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
677  u, v, w, c_xh%B, rho, dt, &
678  ext_bdf%diffusion_coeffs, ext_bdf%ndiff, n)
679  end if
680 
681 
682  call ulag%update()
683  call vlag%update()
684  call wlag%update()
685 
689  call this%user_field_bc_vel%update(this%user_field_bc_vel%field_list, &
690  this%user_field_bc_vel%bc_list, this%c_Xh, t, tstep, "fluid")
691 
692  call this%bc_apply_vel(t, tstep)
693  call this%bc_apply_prs(t, tstep)
694 
695  ! Update material properties if necessary
696  call this%update_material_properties()
697 
698  ! Compute pressure.
699  call profiler_start_region('Pressure_residual', 18)
700  call prs_res%compute(p, p_res,&
701  u, v, w, &
702  u_e, v_e, w_e, &
703  f_x, f_y, f_z, &
704  c_xh, gs_xh, &
705  this%bc_prs_surface, this%bc_sym_surface,&
706  ax_prs, ext_bdf%diffusion_coeffs(1), dt, &
707  mu_field, rho_field)
708 
709  call gs_xh%op(p_res, gs_op_add)
710  call bc_list_apply_scalar(this%bclst_dp, p_res%x, p%dof%size(), t, tstep)
711  call profiler_end_region('Pressure_residual', 18)
712 
713  call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
714  'Pressure')
715 
716  call this%pc_prs%update()
717  call profiler_start_region('Pressure_solve', 3)
718  ksp_results(1) = &
719  this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, this%bclst_dp, gs_xh)
720 
721  call profiler_end_region('Pressure_solve', 3)
722 
723  call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
724  this%bclst_dp, gs_xh, n, tstep, dt_controller)
725 
726  call field_add2(p, dp, n)
727 
728  ! Compute velocity.
729  call profiler_start_region('Velocity_residual', 19)
730  call vel_res%compute(ax_vel, u, v, w, &
731  u_res, v_res, w_res, &
732  p, &
733  f_x, f_y, f_z, &
734  c_xh, msh, xh, &
735  mu_field, rho_field, ext_bdf%diffusion_coeffs(1), &
736  dt, dm_xh%size())
737 
738  call gs_xh%op(u_res, gs_op_add)
739  call gs_xh%op(v_res, gs_op_add)
740  call gs_xh%op(w_res, gs_op_add)
741 
742  call bc_list_apply_vector(this%bclst_vel_res,&
743  u_res%x, v_res%x, w_res%x, dm_xh%size(),&
744  t, tstep)
745 
746  ! We should implement a bc that takes three field_bcs and implements
747  ! vector_apply
748  if (neko_bcknd_device .eq. 1) then
749  call this%bc_field_dirichlet_u%apply_scalar_dev(u_res%x_d, t, tstep)
750  call this%bc_field_dirichlet_v%apply_scalar_dev(v_res%x_d, t, tstep)
751  call this%bc_field_dirichlet_w%apply_scalar_dev(w_res%x_d, t, tstep)
752  else
753  call this%bc_field_dirichlet_u%apply_scalar(u_res%x, n, t, tstep)
754  call this%bc_field_dirichlet_v%apply_scalar(v_res%x, n, t, tstep)
755  call this%bc_field_dirichlet_w%apply_scalar(w_res%x, n, t, tstep)
756  end if
757 
758  call profiler_end_region('Velocity_residual', 19)
759 
760  call this%proj_u%pre_solving(u_res%x, tstep, c_xh, n, dt_controller)
761  call this%proj_v%pre_solving(v_res%x, tstep, c_xh, n, dt_controller)
762  call this%proj_w%pre_solving(w_res%x, tstep, c_xh, n, dt_controller)
763 
764  call this%pc_vel%update()
765 
766  call profiler_start_region("Velocity_solve", 4)
767  ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
768  u_res%x, v_res%x, w_res%x, n, c_xh, &
769  this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh)
770  call profiler_end_region("Velocity_solve", 4)
771 
772  call this%proj_u%post_solving(du%x, ax_vel, c_xh, &
773  this%bclst_du, gs_xh, n, tstep, dt_controller)
774  call this%proj_v%post_solving(dv%x, ax_vel, c_xh, &
775  this%bclst_dv, gs_xh, n, tstep, dt_controller)
776  call this%proj_w%post_solving(dw%x, ax_vel, c_xh, &
777  this%bclst_dw, gs_xh, n, tstep, dt_controller)
778 
779  if (neko_bcknd_device .eq. 1) then
780  call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
781  du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
782  else
783  call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
784  end if
785 
786  if (this%forced_flow_rate) then
787  call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
788  c_xh, gs_xh, ext_bdf, rho, mu,&
789  dt, this%bclst_dp, this%bclst_du, this%bclst_dv, &
790  this%bclst_dw, this%bclst_vel_res, ax_vel, this%ksp_prs, &
791  this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
792  this%ksp_vel%max_iter)
793  end if
794 
795  call fluid_step_info(tstep, t, dt, ksp_results)
796 
797  call this%scratch%relinquish_field(temp_indices)
798 
799  end associate
800  call profiler_end_region('Fluid', 1)
801  end subroutine fluid_pnpn_step
802 
803 
804 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:53
Retrieves a parameter by name or throws an error.
Definition: json_utils.f90:44
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:501
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:583
subroutine, public bc_list_init(bclst, size)
Constructor for a list of boundary conditions.
Definition: bc.f90:460
subroutine, public bc_list_free(bclst)
Destructor for a list of boundary conditions.
Definition: bc.f90:487
subroutine, public bc_list_apply_scalar(bclst, x, n, t, tstep)
Apply a list of boundary conditions to a scalar field.
Definition: bc.f90:526
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
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_step(this, t, tstep, dt, ext_bdf, dt_controller)
Advance fluid simulation in time.
Definition: fluid_pnpn.f90:585
subroutine fluid_pnpn_free(this)
Definition: fluid_pnpn.f90:506
subroutine fluid_pnpn_init(this, msh, lx, params, user, material_properties, time_scheme)
Definition: fluid_pnpn.f90:147
subroutine fluid_pnpn_restart(this, dtlag, tlag)
Definition: fluid_pnpn.f90:369
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
Implements material_properties_t type.
Definition: math.f90:60
subroutine, public col2(a, b, n)
Vector multiplication .
Definition: math.f90:689
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
Contains all the material properties necessary in the simulation.
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.