Neko 1.99.3
A portable framework for high-order spectral element flow simulations
Loading...
Searching...
No Matches
fluid_pnpn.f90
Go to the documentation of this file.
1! Copyright (c) 2022-2025, 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!
35 use, intrinsic :: iso_fortran_env, only : error_unit
36 use coefs, only : coef_t
37 use symmetry, only : symmetry_t
38 use registry, only : neko_registry
39 use logger, only : neko_log, log_size
40 use num_types, only : rp
41 use krylov, only : ksp_monitor_t
43 pnpn_prs_res_factory, pnpn_vel_res_factory, &
44 pnpn_prs_res_stress_factory, pnpn_vel_res_stress_factory
46 rhs_maker_oifs_t, rhs_maker_sumab_fctry, rhs_maker_bdf_fctry, &
47 rhs_maker_ext_fctry, rhs_maker_oifs_fctry
52 use fluid_aux, only : fluid_step_info
53 use projection, only : projection_t
57 use advection, only : advection_t, advection_factory
59 use json_module, only : json_file, json_core, json_value
62 use ax_product, only : ax_t, ax_helm_factory
63 use field, only : field_t
64 use dirichlet, only : dirichlet_t
68 use non_normal, only : non_normal_t
69 use checkpoint, only : chkp_t
70 use mesh, only : mesh_t
71 use user_intf, only : user_t
73 use gs_ops, only : gs_op_add
75 use mathops, only : opadd2cm, opcolv
76 use bc_list, only : bc_list_t
80 use bc, only : bc_t
81 use file, only : file_t
82 use operators, only : ortho, rotate_cyc
83 use opr_device, only : device_ortho
84 use time_state, only : time_state_t
85 use comm, only : neko_comm
86 use ale_manager, only : ale_manager_t
88 use mpi_f08, only : mpi_allreduce, mpi_in_place, mpi_max, mpi_lor, &
89 mpi_integer, mpi_logical
90 implicit none
91 private
92
93
94
96
99 integer :: schwarz_iterations = 0
100
102 type(field_t) :: p_res, u_res, v_res, w_res
103
106 type(field_t) :: dp, du, dv, dw
107
109 type(ale_manager_t) :: ale
110
111 ! ! Implicit operators, i.e. the left-hand-side of the Helmholz problem.
112 !
113
114 ! Coupled Helmholz operator for velocity
115 class(ax_t), allocatable :: ax_vel
116 ! Helmholz operator for pressure
117 class(ax_t), allocatable :: ax_prs
118
119 !
120 ! Projections for solver speed-up
121 !
122
124 type(projection_t) :: proj_prs
125 type(projection_vel_t) :: proj_vel
126
127 !
128 ! Special Karniadakis scheme boundary conditions in the pressure equation
129 !
130
132 type(facet_normal_t) :: bc_prs_surface
133
135 type(facet_normal_t) :: bc_sym_surface
136
137 !
138 ! Boundary conditions and lists for residuals and solution increments
139 !
140
142 type(zero_dirichlet_t) :: bc_vel_res
144 type(zero_dirichlet_t) :: bc_du
146 type(zero_dirichlet_t) :: bc_dv
148 type(zero_dirichlet_t) :: bc_dw
150 type(zero_dirichlet_t) :: bc_dp
151
153 type(bc_list_t) :: bclst_vel_res
154 type(bc_list_t) :: bclst_du
155 type(bc_list_t) :: bclst_dv
156 type(bc_list_t) :: bclst_dw
157 type(bc_list_t) :: bclst_dp
158
159
160 ! Checker for wether we have a strong pressure bc. If not, the pressure
161 ! is demeaned at every time step.
162 logical :: prs_dirichlet = .false.
163
164
165 ! The advection operator.
166 class(advection_t), allocatable :: adv
167
168 ! Time OIFS interpolation scheme for advection.
169 logical :: oifs
170
171 ! Time variables
172 type(field_t) :: abx1, aby1, abz1
173 type(field_t) :: abx2, aby2, abz2
174
175 ! Advection terms for the oifs method
176 type(field_t) :: advx, advy, advz
177
179 class(pnpn_prs_res_t), allocatable :: prs_res
180
182 class(pnpn_vel_res_t), allocatable :: vel_res
183
185 class(rhs_maker_sumab_t), allocatable :: sumab
186
188 class(rhs_maker_ext_t), allocatable :: makeabf
189
191 class(rhs_maker_bdf_t), allocatable :: makebdf
192
194 class(rhs_maker_oifs_t), allocatable :: makeoifs
195
197 type(fluid_volflow_t) :: vol_flow
198
200 logical :: full_stress_formulation = .false.
201
202 contains
204 procedure, pass(this) :: init => fluid_pnpn_init
206 procedure, pass(this) :: free => fluid_pnpn_free
208 procedure, pass(this) :: step => fluid_pnpn_step
210 procedure, pass(this) :: restart => fluid_pnpn_restart
212 procedure, pass(this) :: setup_bcs => fluid_pnpn_setup_bcs
214 procedure, pass(this) :: write_boundary_conditions => &
216 end type fluid_pnpn_t
217
218 interface
219
226 module subroutine pressure_bc_factory(object, scheme, json, coef, user)
227 class(bc_t), pointer, intent(inout) :: object
228 type(fluid_pnpn_t), intent(in) :: scheme
229 type(json_file), intent(inout) :: json
230 type(coef_t), target, intent(in) :: coef
231 type(user_t), intent(in) :: user
232 end subroutine pressure_bc_factory
233 end interface
234
235 interface
236
243 module subroutine velocity_bc_factory(object, scheme, json, coef, user)
244 class(bc_t), pointer, intent(inout) :: object
245 type(fluid_pnpn_t), intent(inout) :: scheme
246 type(json_file), intent(inout) :: json
247 type(coef_t), target, intent(in) :: coef
248 type(user_t), intent(in) :: user
249 end subroutine velocity_bc_factory
250 end interface
251
252contains
253
254 subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
255 class(fluid_pnpn_t), target, intent(inout) :: this
256 type(mesh_t), target, intent(inout) :: msh
257 integer, intent(in) :: lx
258 type(json_file), target, intent(inout) :: params
259 type(user_t), target, intent(in) :: user
260 type(chkp_t), target, intent(inout) :: chkp
261 character(len=15), parameter :: scheme = 'Modular (Pn/Pn)'
262 integer :: i
263 class(bc_t), pointer :: bc_i, vel_bc
264 real(kind=rp) :: abs_tol
265 character(len=LOG_SIZE) :: log_buf
266 integer :: ierr, integer_val, solver_maxiter
267 character(len=:), allocatable :: solver_type, precon_type
268 logical :: monitor, found
269 logical :: advection
270 type(json_file) :: numerics_params, precon_params
271
272 call this%free()
273
274 ! Initialize base class
275 call this%init_base(msh, lx, params, scheme, user, .true.)
276
277 ! Add pressure field to the registry. For this scheme it is in the same
278 ! Xh as the velocity
279 call neko_registry%add_field(this%dm_Xh, 'p')
280 this%p => neko_registry%get_field('p')
281
282 !
283 ! Select governing equations via associated residual and Ax types
284 !
285
286 call json_get_or_lookup(params, 'case.numerics.time_order', integer_val)
287 allocate(this%ext_bdf)
288 call this%ext_bdf%init(integer_val)
289
290 call json_get_or_default(params, "case.fluid.full_stress_formulation", &
291 this%full_stress_formulation, .false.)
292
293 call json_get_or_default(params, "case.fluid.cyclic", this%c_Xh%cyclic, &
294 .false.)
295 call this%c_Xh%generate_cyclic_bc()
296
297 if (this%full_stress_formulation) then
298 ! Setup backend dependent Ax routines
299 call ax_helm_factory(this%Ax_vel, full_formulation = .true.)
300
301 ! Setup backend dependent prs residual routines
302 call pnpn_prs_res_stress_factory(this%prs_res)
303
304 ! Setup backend dependent vel residual routines
305 call pnpn_vel_res_stress_factory(this%vel_res)
306 else
307 ! Setup backend dependent Ax routines
308 call ax_helm_factory(this%Ax_vel, full_formulation = .false.)
309
310 ! Setup backend dependent prs residual routines
311 call pnpn_prs_res_factory(this%prs_res)
312
313 ! Setup backend dependent vel residual routines
314 call pnpn_vel_res_factory(this%vel_res)
315 end if
316
317
318
319 if (params%valid_path('case.fluid.nut_field')) then
320 if (.not. this%full_stress_formulation) then
321 call neko_error("You need to set full_stress_formulation to " // &
322 "true for the fluid to have a spatially varying " // &
323 "viscocity field.")
324 end if
325 call json_get(params, 'case.fluid.nut_field', this%nut_field_name)
326 else
327 this%nut_field_name = ""
328 end if
329
330 ! Setup Ax for the pressure
331 call ax_helm_factory(this%Ax_prs, full_formulation = .false.)
332
333
334 ! Setup backend dependent summation of AB/BDF
335 call rhs_maker_sumab_fctry(this%sumab)
336
337 ! Setup backend dependent summation of extrapolation scheme
338 call rhs_maker_ext_fctry(this%makeabf)
339
340 ! Setup backend depenent contributions to F from lagged BD terms
341 call rhs_maker_bdf_fctry(this%makebdf)
342
343 ! Setup backend dependent summations of the OIFS method
344 call rhs_maker_oifs_fctry(this%makeoifs)
345
346 ! Initialize variables specific to this plan
347 associate(xh_lx => this%Xh%lx, xh_ly => this%Xh%ly, xh_lz => this%Xh%lz, &
348 dm_xh => this%dm_Xh, nelv => this%msh%nelv)
349
350 call this%p_res%init(dm_xh, "p_res")
351 call this%u_res%init(dm_xh, "u_res")
352 call this%v_res%init(dm_xh, "v_res")
353 call this%w_res%init(dm_xh, "w_res")
354 call this%abx1%init(dm_xh, "abx1")
355 call this%aby1%init(dm_xh, "aby1")
356 call this%abz1%init(dm_xh, "abz1")
357 call this%abx2%init(dm_xh, "abx2")
358 call this%aby2%init(dm_xh, "aby2")
359 call this%abz2%init(dm_xh, "abz2")
360 call this%advx%init(dm_xh, "advx")
361 call this%advy%init(dm_xh, "advy")
362 call this%advz%init(dm_xh, "advz")
363 end associate
364
365 call this%du%init(this%dm_Xh, 'du')
366 call this%dv%init(this%dm_Xh, 'dv')
367 call this%dw%init(this%dm_Xh, 'dw')
368 call this%dp%init(this%dm_Xh, 'dp')
369 ! Initialize ALE
370 call this%ale%init(this%c_Xh, params, user)
371
372 call neko_log%section("Fluid boundary conditions")
373 ! Set up boundary conditions
374 call this%setup_bcs(user, params)
375
376 ! Check if we need to output boundaries
377 call json_get_or_default(params, 'case.output_boundary', found, .false.)
378 if (found) call this%write_boundary_conditions()
379 call neko_log%end_section()
380
381 call this%proj_prs%init(this%dm_Xh%size(), this%pr_projection_dim, &
382 this%pr_projection_activ_step, &
383 this%pr_projection_reorthogonalize_basis)
384
385 call this%proj_vel%init(this%dm_Xh%size(), this%vel_projection_dim, &
386 this%vel_projection_activ_step)
387
388
389
390 ! Determine the time-interpolation scheme
391 call json_get_or_default(params, 'case.numerics.oifs', this%oifs, .false.)
392 if (params%valid_path('case.fluid.flow_rate_force')) then
393 call this%vol_flow%init(this%dm_Xh, params)
394 end if
395
396 ! Setup pressure solver
397 call neko_log%section("Pressure solver")
398
399 call json_get_or_lookup_or_default(params, &
400 'case.fluid.pressure_solver.max_iterations', &
401 solver_maxiter, 800)
402 call json_get(params, 'case.fluid.pressure_solver.type', solver_type)
403 call json_get(params, 'case.fluid.pressure_solver.preconditioner.type', &
404 precon_type)
405 call json_get(params, &
406 'case.fluid.pressure_solver.preconditioner', precon_params)
407 call json_get_or_lookup(params, &
408 'case.fluid.pressure_solver.absolute_tolerance', &
409 abs_tol)
410 call json_get_or_default(params, 'case.fluid.pressure_solver.monitor', &
411 monitor, .false.)
412 call neko_log%message('Type : ('// trim(solver_type) // &
413 ', ' // trim(precon_type) // ')')
414 write(log_buf, '(A,ES13.6)') 'Abs tol :', abs_tol
415 call neko_log%message(log_buf)
416
417 call this%solver_factory(this%ksp_prs, this%dm_Xh%size(), &
418 solver_type, solver_maxiter, abs_tol, monitor)
419 call this%precon_factory_(this%pc_prs, this%ksp_prs, &
420 this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_prs, &
421 precon_type, precon_params)
422 call neko_log%end_section()
423
424 ! Initialize the advection factory
425 call json_get_or_default(params, 'case.fluid.advection', advection, .true.)
426 call json_get(params, 'case.numerics', numerics_params)
427 call advection_factory(this%adv, numerics_params, this%c_Xh, &
428 this%ulag, this%vlag, this%wlag, &
429 chkp%dtlag, chkp%tlag, this%ext_bdf, &
430 .not. advection)
431 ! Should be in init_base maybe?
432 this%chkp => chkp
433 ! This is probably scheme specific
434 call this%chkp%add_fluid(this%u, this%v, this%w, this%p)
435
436 this%chkp%abx1 => this%abx1
437 this%chkp%abx2 => this%abx2
438 this%chkp%aby1 => this%aby1
439 this%chkp%aby2 => this%aby2
440 this%chkp%abz1 => this%abz1
441 this%chkp%abz2 => this%abz2
442 call this%chkp%add_lag(this%ulag, this%vlag, this%wlag)
443
444 ! Add checkpoint data for ALE.
445 if (this%ale%active) then
446 call this%chkp%add_ale(this%c_Xh%dof%x, this%c_Xh%dof%y, &
447 this%c_Xh%dof%z, &
448 this%c_Xh%Blag, this%c_Xh%Blaglag, &
449 this%ale%wm_x, this%ale%wm_y, this%ale%wm_z, &
450 this%ale%wm_x_lag, this%ale%wm_y_lag, &
451 this%ale%wm_z_lag, &
452 this%ale%global_pivot_pos, &
453 this%ale%global_pivot_vel_lag, &
454 this%ale%global_basis_pos, &
455 this%ale%global_basis_vel_lag)
456 end if
457
459 call json_get_or_default(params, 'case.fluid.schwarz_iterations', &
460 this%schwarz_iterations, 0)
461
462 call neko_log%end_section()
463
464 end subroutine fluid_pnpn_init
465
466 subroutine fluid_pnpn_restart(this, chkp)
467 class(fluid_pnpn_t), target, intent(inout) :: this
468 type(chkp_t), intent(inout) :: chkp
469 real(kind=rp) :: dtlag(10), tlag(10)
470 type(field_t) :: u_temp, v_temp, w_temp
471 integer :: i, j, n
472
473 dtlag = chkp%dtlag
474 tlag = chkp%tlag
475
476 n = this%u%dof%size()
477 if (allocated(chkp%previous_mesh%elements) .or. &
478 chkp%previous_Xh%lx .ne. this%Xh%lx) then
479 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
480 c_xh => this%c_Xh, ulag => this%ulag, vlag => this%vlag, &
481 wlag => this%wlag)
482 do concurrent(j = 1:n)
483 u%x(j,1,1,1) = u%x(j,1,1,1) * c_xh%mult(j,1,1,1)
484 v%x(j,1,1,1) = v%x(j,1,1,1) * c_xh%mult(j,1,1,1)
485 w%x(j,1,1,1) = w%x(j,1,1,1) * c_xh%mult(j,1,1,1)
486 p%x(j,1,1,1) = p%x(j,1,1,1) * c_xh%mult(j,1,1,1)
487 end do
488 do i = 1, this%ulag%size()
489 do concurrent(j = 1:n)
490 ulag%lf(i)%x(j,1,1,1) = ulag%lf(i)%x(j,1,1,1) &
491 * c_xh%mult(j,1,1,1)
492 vlag%lf(i)%x(j,1,1,1) = vlag%lf(i)%x(j,1,1,1) &
493 * c_xh%mult(j,1,1,1)
494 wlag%lf(i)%x(j,1,1,1) = wlag%lf(i)%x(j,1,1,1) &
495 * c_xh%mult(j,1,1,1)
496 end do
497 end do
498 end associate
499 end if
500
501 if (neko_bcknd_device .eq. 1) then
502 associate(u => this%u, v => this%v, w => this%w, &
503 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag,&
504 p => this%p)
505 call device_memcpy(u%x, u%x_d, u%dof%size(), &
506 host_to_device, sync = .false.)
507 call device_memcpy(v%x, v%x_d, v%dof%size(), &
508 host_to_device, sync = .false.)
509 call device_memcpy(w%x, w%x_d, w%dof%size(), &
510 host_to_device, sync = .false.)
511 call device_memcpy(p%x, p%x_d, p%dof%size(), &
512 host_to_device, sync = .false.)
513 call device_memcpy(ulag%lf(1)%x, ulag%lf(1)%x_d, &
514 u%dof%size(), host_to_device, sync = .false.)
515 call device_memcpy(ulag%lf(2)%x, ulag%lf(2)%x_d, &
516 u%dof%size(), host_to_device, sync = .false.)
517
518 call device_memcpy(vlag%lf(1)%x, vlag%lf(1)%x_d, &
519 v%dof%size(), host_to_device, sync = .false.)
520 call device_memcpy(vlag%lf(2)%x, vlag%lf(2)%x_d, &
521 v%dof%size(), host_to_device, sync = .false.)
522
523 call device_memcpy(wlag%lf(1)%x, wlag%lf(1)%x_d, &
524 w%dof%size(), host_to_device, sync = .false.)
525 call device_memcpy(wlag%lf(2)%x, wlag%lf(2)%x_d, &
526 w%dof%size(), host_to_device, sync = .false.)
527 call device_memcpy(this%abx1%x, this%abx1%x_d, &
528 w%dof%size(), host_to_device, sync = .false.)
529 call device_memcpy(this%abx2%x, this%abx2%x_d, &
530 w%dof%size(), host_to_device, sync = .false.)
531 call device_memcpy(this%aby1%x, this%aby1%x_d, &
532 w%dof%size(), host_to_device, sync = .false.)
533 call device_memcpy(this%aby2%x, this%aby2%x_d, &
534 w%dof%size(), host_to_device, sync = .false.)
535 call device_memcpy(this%abz1%x, this%abz1%x_d, &
536 w%dof%size(), host_to_device, sync = .false.)
537 call device_memcpy(this%abz2%x, this%abz2%x_d, &
538 w%dof%size(), host_to_device, sync = .false.)
539 call device_memcpy(this%advx%x, this%advx%x_d, &
540 w%dof%size(), host_to_device, sync = .false.)
541 call device_memcpy(this%advy%x, this%advy%x_d, &
542 w%dof%size(), host_to_device, sync = .false.)
543 call device_memcpy(this%advz%x, this%advz%x_d, &
544 w%dof%size(), host_to_device, sync = .false.)
545 end associate
546 end if
547 ! Make sure that continuity is maintained (important for interpolation)
548 ! Do not do this for lagged rhs
549 ! (derivatives are not necessairly coninous across elements)
550
551 if (allocated(chkp%previous_mesh%elements) &
552 .or. chkp%previous_Xh%lx .ne. this%Xh%lx) then
553
554 call rotate_cyc(this%u, this%v, this%w, 1, this%c_Xh)
555 call this%gs_Xh%op(this%u, gs_op_add)
556 call this%gs_Xh%op(this%v, gs_op_add)
557 call this%gs_Xh%op(this%w, gs_op_add)
558 call this%gs_Xh%op(this%p, gs_op_add)
559 call rotate_cyc(this%u, this%v, this%w, 0, this%c_Xh)
560
561 do i = 1, this%ulag%size()
562 call rotate_cyc(this%ulag%lf(i), this%vlag%lf(i), &
563 this%wlag%lf(i), 1, this%c_Xh)
564 call this%gs_Xh%op(this%ulag%lf(i), gs_op_add)
565 call this%gs_Xh%op(this%vlag%lf(i), gs_op_add)
566 call this%gs_Xh%op(this%wlag%lf(i), gs_op_add)
567 call rotate_cyc(this%ulag%lf(i), this%vlag%lf(i), &
568 this%wlag%lf(i), 0, this%c_Xh)
569 end do
570 end if
571
572 call this%ale%set_coef_restart(this%c_Xh, this%adv, chkp%t)
573
574
575 end subroutine fluid_pnpn_restart
576
577 subroutine fluid_pnpn_free(this)
578 class(fluid_pnpn_t), intent(inout) :: this
579
580 !Deallocate velocity and pressure fields
581 call this%scheme_free()
582
583 if (allocated(this%ext_bdf)) then
584 call this%ext_bdf%free()
585 deallocate(this%ext_bdf)
586 end if
587
588 call this%bc_vel_res%free()
589 call this%bc_du%free()
590 call this%bc_dv%free()
591 call this%bc_dw%free()
592 call this%bc_dp%free()
593
594 call this%bc_prs_surface%free()
595 call this%bc_sym_surface%free()
596 call this%bclst_vel_res%free()
597 call this%bclst_du%free()
598 call this%bclst_dv%free()
599 call this%bclst_dw%free()
600 call this%bclst_dp%free()
601 call this%proj_prs%free()
602 call this%proj_vel%free()
603
604 call this%p_res%free()
605 call this%u_res%free()
606 call this%v_res%free()
607 call this%w_res%free()
608
609 call this%ale%free()
610
611 call this%du%free()
612 call this%dv%free()
613 call this%dw%free()
614 call this%dp%free()
615
616 call this%abx1%free()
617 call this%aby1%free()
618 call this%abz1%free()
619
620 call this%abx2%free()
621 call this%aby2%free()
622 call this%abz2%free()
623
624 call this%advx%free()
625 call this%advy%free()
626 call this%advz%free()
627
628 if (allocated(this%adv)) then
629 call this%adv%free()
630 deallocate(this%adv)
631 end if
632
633 if (allocated(this%Ax_vel)) then
634 deallocate(this%Ax_vel)
635 end if
636
637 if (allocated(this%Ax_prs)) then
638 deallocate(this%Ax_prs)
639 end if
640
641 if (allocated(this%prs_res)) then
642 deallocate(this%prs_res)
643 end if
644
645 if (allocated(this%vel_res)) then
646 deallocate(this%vel_res)
647 end if
648
649 if (allocated(this%sumab)) then
650 deallocate(this%sumab)
651 end if
652
653 if (allocated(this%makeabf)) then
654 deallocate(this%makeabf)
655 end if
656
657 if (allocated(this%makebdf)) then
658 deallocate(this%makebdf)
659 end if
660
661 if (allocated(this%makeoifs)) then
662 deallocate(this%makeoifs)
663 end if
664
665 if (allocated(this%ext_bdf)) then
666 deallocate(this%ext_bdf)
667 end if
668
669 call this%vol_flow%free()
670
671 end subroutine fluid_pnpn_free
672
679 subroutine fluid_pnpn_step(this, time, dt_controller)
680 class(fluid_pnpn_t), target, intent(inout) :: this
681 type(time_state_t), intent(in) :: time
682 type(time_step_controller_t), intent(in) :: dt_controller
683 ! number of degrees of freedom
684 integer :: n
685 ! Solver results monitors (pressure + 3 velocity)
686 type(ksp_monitor_t) :: ksp_results(4)
687 integer :: iter
688
689 type(file_t) :: dump_file
690 class(bc_t), pointer :: bc_i
691 type(non_normal_t), pointer :: bc_j
692
693 if (this%freeze) return
694
695 n = this%dm_Xh%size()
696
697 call profiler_start_region('Fluid', 1)
698 associate(u => this%u, v => this%v, w => this%w, p => this%p, &
699 u_e => this%u_e, v_e => this%v_e, w_e => this%w_e, &
700 du => this%du, dv => this%dv, dw => this%dw, dp => this%dp, &
701 u_res => this%u_res, v_res => this%v_res, w_res => this%w_res, &
702 p_res => this%p_res, ax_vel => this%Ax_vel, ax_prs => this%Ax_prs, &
703 xh => this%Xh, &
704 c_xh => this%c_Xh, dm_xh => this%dm_Xh, gs_xh => this%gs_Xh, &
705 ulag => this%ulag, vlag => this%vlag, wlag => this%wlag, &
706 msh => this%msh, prs_res => this%prs_res, &
707 source_term => this%source_term, vel_res => this%vel_res, &
708 sumab => this%sumab, makeoifs => this%makeoifs, &
709 makeabf => this%makeabf, makebdf => this%makebdf, &
710 vel_projection_dim => this%vel_projection_dim, &
711 pr_projection_dim => this%pr_projection_dim, &
712 oifs => this%oifs, &
713 rho => this%rho, mu_tot => this%mu_tot, &
714 f_x => this%f_x, f_y => this%f_y, f_z => this%f_z, &
715 t => time%t, tstep => time%tstep, dt => time%dt, &
716 ext_bdf => this%ext_bdf, event => glb_cmd_event, &
717 ale => this%ale)
718
719 ! Extrapolate the velocity if it's not done in nut_field estimation
720 call sumab%compute_fluid(u_e, v_e, w_e, u, v, w, &
721 ulag, vlag, wlag, ext_bdf%advection_coeffs%x, ext_bdf%nadv)
722
723 ! Compute the source terms
724 call this%source_term%compute(time)
725
726 ! Add Neumann bc contributions to the RHS
727 call this%bcs_vel%apply_vector(f_x%x, f_y%x, f_z%x, &
728 this%dm_Xh%size(), time, strong = .false.)
729
730 if (this%ale%active) then
731 if (oifs) then
732 call neko_error("ALE is not yet supported " // &
733 "with OIFS time integration.")
734 end if
736 call this%adv%compute_ale(u, v, w, &
737 ale%wm_x, ale%wm_y, ale%wm_z, &
738 f_x, f_y, f_z, &
739 xh, c_xh, dm_xh%size())
740 end if
741
742
743 if (oifs) then
744 ! Add the advection operators to the right-hand-side.
745 call this%adv%compute(u, v, w, &
746 this%advx, this%advy, this%advz, &
747 xh, this%c_Xh, dm_xh%size(), dt)
748
749 ! At this point the RHS contains the sum of the advection operator and
750 ! additional source terms, evaluated using the velocity field from the
751 ! previous time-step. Now, this value is used in the explicit time
752 ! scheme to advance both terms in time.
753
754 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
755 this%abx2, this%aby2, this%abz2, &
756 f_x%x, f_y%x, f_z%x, &
757 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
758
759 ! Now, the source terms from the previous time step are added to the RHS.
760 call makeoifs%compute_fluid(this%advx%x, this%advy%x, this%advz%x, &
761 f_x%x, f_y%x, f_z%x, &
762 rho%x(1,1,1,1), dt, n)
763 else
764 ! Add the advection operators to the right-hand-side.
765 call this%adv%compute(u, v, w, &
766 f_x, f_y, f_z, &
767 xh, this%c_Xh, dm_xh%size())
768
769 ! At this point the RHS contains the sum of the advection operator and
770 ! additional source terms, evaluated using the velocity field from the
771 ! previous time-step. Now, this value is used in the explicit time
772 ! scheme to advance both terms in time.
773
774 call makeabf%compute_fluid(this%abx1, this%aby1, this%abz1,&
775 this%abx2, this%aby2, this%abz2, &
776 f_x%x, f_y%x, f_z%x, &
777 rho%x(1,1,1,1), ext_bdf%advection_coeffs%x, n)
778
779 ! Add the RHS contributions coming from the BDF scheme.
780 ! Blag and Blaglag are history of B matrices, mainly used for ALE.
781 ! For a normal simulation (no moving mesh), Blag and Blaglag
782 ! are just the initial B matrix, filled at initialization.
783 call makebdf%compute_fluid(ulag, vlag, wlag, f_x%x, f_y%x, f_z%x, &
784 u, v, w, c_xh%B, rho%x(1,1,1,1), dt, &
785 ext_bdf%diffusion_coeffs%x, ext_bdf%ndiff, n, &
786 c_xh%Blag, c_xh%Blaglag)
787
788 end if
789
790 if (this%ale%active) then
791 ! Advance Mesh (Moves points, updates B history, updates wm_lags)
792 call this%ale%advance_mesh(c_xh, time, ext_bdf%nadv)
793
794 call profiler_start_region('ALE recompute metrics')
795 ! Update Metrics
796 call c_xh%recompute_metrics()
797 ! Update the metrics used by the adv operator for delaiasing (coef_GL)
798 ! Maps the updated coef_GLL to coef_GL.
799 call this%adv%recompute_metrics(c_xh, .true.)
800 call profiler_end_region('ALE recompute metrics')
801 end if
802
803 call ulag%update()
804 call vlag%update()
805 call wlag%update()
806
807 ! Update material properties if necessary
808 call this%update_material_properties(time)
809
810 do iter = 1, 1 + this%schwarz_iterations
811
812 call this%bc_apply_vel(time, strong = .true.)
813 call this%bc_apply_prs(time)
814
815 ! Compute pressure residual.
816 call profiler_start_region('Pressure_residual', 18)
817 call prs_res%compute(p, p_res,&
818 u, v, w, &
819 u_e, v_e, w_e, &
820 f_x, f_y, f_z, &
821 c_xh, gs_xh, &
822 this%bc_prs_surface, this%bc_sym_surface,&
823 ax_prs, ext_bdf%diffusion_coeffs%x(1), dt, &
824 mu_tot, rho, event)
825
826
827 ! De-mean the pressure residual when no strong pressure boundaries present
828 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1) then
829 call device_ortho(p_res%x_d, this%glb_n_points, n)
830 else if (.not. this%prs_dirichlet) then
831 call ortho(p_res%x, this%glb_n_points, n)
832 end if
833
834 call gs_xh%op(p_res, gs_op_add, event)
835 call device_event_sync(event)
836
837 ! Set the residual to zero at strong pressure boundaries.
838 call this%bclst_dp%apply_scalar(p_res%x, p%dof%size(), time)
839
840
841 call profiler_end_region('Pressure_residual', 18)
842
843 ! Do projections only on the actual solutions of the tstep
844 ! not intermediate solutions from the subiterations.
845 if (iter .eq. 1) then
846 call this%proj_prs%pre_solving(p_res%x, tstep, c_xh, n, dt_controller, &
847 ax = ax_prs, gs_h = gs_xh, bclst = this%bclst_dp, &
848 string = 'Pressure')
849 end if
850
851 call this%pc_prs%update()
852
853 call profiler_start_region('Pressure_solve', 3)
854
855 ! Solve for the pressure increment.
856 ksp_results(1) = &
857 this%ksp_prs%solve(ax_prs, dp, p_res%x, n, c_xh, &
858 this%bclst_dp, gs_xh)
859 ksp_results(1)%name = 'Pressure'
860
861
862 call profiler_end_region('Pressure_solve', 3)
863
864 if (iter .eq. 1) then
865 call this%proj_prs%post_solving(dp%x, ax_prs, c_xh, &
866 this%bclst_dp, gs_xh, n, tstep, dt_controller)
867 end if
868
869 ! Update the pressure with the increment. Demean if necessary.
870 call field_add2(p, dp, n)
871 if (.not. this%prs_dirichlet .and. neko_bcknd_device .eq. 1) then
872 call device_ortho(p%x_d, this%glb_n_points, n)
873 else if (.not. this%prs_dirichlet) then
874 call ortho(p%x, this%glb_n_points, n)
875 end if
876
877 ! Compute velocity residual.
878 call profiler_start_region('Velocity_residual', 19)
879 call vel_res%compute(ax_vel, u, v, w, &
880 u_res, v_res, w_res, &
881 p, &
882 f_x, f_y, f_z, &
883 c_xh, msh, xh, &
884 mu_tot, rho, ext_bdf%diffusion_coeffs%x(1), &
885 dt, dm_xh%size())
886
887 call rotate_cyc(u_res, v_res, w_res, 1, c_xh)
888 call gs_xh%op(u_res, gs_op_add, event)
889 call device_event_sync(event)
890 call gs_xh%op(v_res, gs_op_add, event)
891 call device_event_sync(event)
892 call gs_xh%op(w_res, gs_op_add, event)
893 call device_event_sync(event)
894 call rotate_cyc(u_res, v_res, w_res, 0, c_xh)
895
896 ! Set residual to zero at strong velocity boundaries.
897 call this%bclst_vel_res%apply(u_res, v_res, w_res, time)
898
899
900 call profiler_end_region('Velocity_residual', 19)
901
902 if (iter .eq. 1) then
903 call this%proj_vel%pre_solving(u_res%x, v_res%x, w_res%x, &
904 tstep, c_xh, n, dt_controller, 'Velocity')
905 end if
906
907 call this%pc_vel%update()
908
909 call profiler_start_region("Velocity_solve", 4)
910 ksp_results(2:4) = this%ksp_vel%solve_coupled(ax_vel, du, dv, dw, &
911 u_res%x, v_res%x, w_res%x, n, c_xh, &
912 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, &
913 this%ksp_vel%max_iter)
914 call profiler_end_region("Velocity_solve", 4)
915 if (this%full_stress_formulation) then
916 ksp_results(2)%name = 'Momentum'
917 else
918 ksp_results(2)%name = 'X-Velocity'
919 ksp_results(3)%name = 'Y-Velocity'
920 ksp_results(4)%name = 'Z-Velocity'
921 end if
922
923 if (iter .eq. 1) then
924 call this%proj_vel%post_solving(du%x, dv%x, dw%x, ax_vel, c_xh, &
925 this%bclst_du, this%bclst_dv, this%bclst_dw, gs_xh, n, tstep, &
926 dt_controller)
927 end if
928
929 if (neko_bcknd_device .eq. 1) then
930 call device_opadd2cm(u%x_d, v%x_d, w%x_d, &
931 du%x_d, dv%x_d, dw%x_d, 1.0_rp, n, msh%gdim)
932 else
933 call opadd2cm(u%x, v%x, w%x, du%x, dv%x, dw%x, 1.0_rp, n, msh%gdim)
934 end if
935
936 call fluid_step_info(time, ksp_results, &
937 this%full_stress_formulation, this%strict_convergence, &
938 this%allow_stabilization, iter)
939
940 end do
941
942 if (this%forced_flow_rate) then
943 ! Horrible mu hack?!
944 call this%vol_flow%adjust( u, v, w, p, u_res, v_res, w_res, p_res, &
945 c_xh, gs_xh, ext_bdf, rho%x(1,1,1,1), mu_tot, &
946 dt, time, this%bclst_dp, this%bclst_du, this%bclst_dv, &
947 this%bclst_dw, this%bclst_vel_res, ax_vel, ax_prs, this%ksp_prs, &
948 this%ksp_vel, this%pc_prs, this%pc_vel, this%ksp_prs%max_iter, &
949 this%ksp_vel%max_iter)
950 end if
951
952 ! Update mesh velocities for ALE
953 ! We update them here (end of step) for the next step.
954 ! Returns if .not. ale.
955 call this%ale%update_mesh_velocity(c_xh, time)
956
957 end associate
958 call profiler_end_region('Fluid', 1)
959 end subroutine fluid_pnpn_step
960
963 subroutine fluid_pnpn_setup_bcs(this, user, params)
964 class(fluid_pnpn_t), target, intent(inout) :: this
965 type(user_t), target, intent(in) :: user
966 type(json_file), intent(inout) :: params
967 integer :: i, n_bcs, zone_index, j, zone_size, global_zone_size, ierr
968 class(bc_t), pointer :: bc_i
969 type(json_core) :: core
970 type(json_value), pointer :: bc_object
971 type(json_file) :: bc_subdict
972 logical :: ale_active_local, any_moving_wall, moving_
973 logical :: found
974 ! Monitor which boundary zones have been marked
975 logical, allocatable :: marked_zones(:)
976 integer, allocatable :: zone_indices(:)
977
978 ! For ALE, we set a flag while reading the BCs
979 character(len=:), allocatable :: bc_type_str
980 this%ale%has_moving_boundary = .false.
981 any_moving_wall = .false.
982 ale_active_local = .false.
983 call json_get_or_default(params, &
984 'case.fluid.ale.enabled', &
985 ale_active_local, .false.)
986
987 ! Lists for the residuals and solution increments
988 call this%bclst_vel_res%init()
989 call this%bclst_du%init()
990 call this%bclst_dv%init()
991 call this%bclst_dw%init()
992 call this%bclst_dp%init()
993
994 call this%bc_vel_res%init_from_components(this%c_Xh)
995 call this%bc_du%init_from_components(this%c_Xh)
996 call this%bc_dv%init_from_components(this%c_Xh)
997 call this%bc_dw%init_from_components(this%c_Xh)
998 call this%bc_dp%init_from_components(this%c_Xh)
999
1000 ! Special PnPn boundary conditions for pressure
1001 call this%bc_prs_surface%init_from_components(this%c_Xh)
1002 call this%bc_sym_surface%init_from_components(this%c_Xh)
1003
1004 ! Populate bcs_vel and bcs_prs based on the case file
1005 if (params%valid_path('case.fluid.boundary_conditions')) then
1006 call params%info('case.fluid.boundary_conditions', n_children = n_bcs)
1007 call params%get_core(core)
1008 call params%get('case.fluid.boundary_conditions', bc_object, found)
1009
1010 !
1011 ! Velocity bcs
1012 !
1013 call this%bcs_vel%init(n_bcs)
1014
1015 allocate(marked_zones(size(this%msh%labeled_zones)))
1016 marked_zones = .false.
1017
1018 do i = 1, n_bcs
1019 ! Create a new json containing just the subdict for this bc
1020 call json_extract_item(core, bc_object, i, bc_subdict)
1021
1022 call json_get_or_lookup(bc_subdict, "zone_indices", zone_indices)
1023
1024 ! Set the ALE flag to true if there is any moving no_slip wall
1025 call json_get(bc_subdict, "type", bc_type_str)
1026 moving_ = .false.
1027 if (trim(bc_type_str) .eq. "no_slip") then
1028 call json_get_or_default(bc_subdict, "moving", moving_, .false.)
1029 end if
1030 if (moving_) then
1031 this%ale%has_moving_boundary = .true.
1032 end if
1033
1034 ! Check that we are not trying to assing a bc to zone, for which one
1035 ! has already been assigned and that the zone has more than 0 size
1036 ! in the mesh.
1037 do j = 1, size(zone_indices)
1038 zone_size = this%msh%labeled_zones(zone_indices(j))%size
1039 call mpi_allreduce(zone_size, global_zone_size, 1, &
1040 mpi_integer, mpi_max, neko_comm, ierr)
1041
1042 if (global_zone_size .eq. 0) then
1043 write(error_unit, '(A, A, I0, A, A, I0, A)') "*** ERROR ***: ",&
1044 "Zone index ", zone_indices(j), &
1045 " is invalid as this zone has 0 size, meaning it ", &
1046 "is not in the mesh. Check fluid boundary condition ", &
1047 i, "."
1048 error stop
1049 end if
1050
1051 if (marked_zones(zone_indices(j))) then
1052 write(error_unit, '(A, A, I0, A, A, A, A)') "*** ERROR ***: ", &
1053 "Zone with index ", zone_indices(j), &
1054 " has already been assigned a boundary condition. ", &
1055 "Please check your boundary_conditions entry for the ", &
1056 "fluid and make sure that each zone index appears only ", &
1057 "in a single boundary condition."
1058 error stop
1059 else
1060 marked_zones(zone_indices(j)) = .true.
1061 end if
1062 end do
1063
1064 bc_i => null()
1065 call velocity_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1066
1067 ! Not all bcs require an allocation for velocity in particular,
1068 ! so we check.
1069 if (associated(bc_i)) then
1070
1071 ! We need to treat mixed bcs separately because they are by
1072 ! convention marked weak and currently contain nested
1073 ! bcs, some of which are strong.
1074 select type (bc_i)
1075 type is (symmetry_t)
1076 ! Symmetry has 3 internal bcs, but only one actually contains
1077 ! markings.
1078 ! Symmetry's apply_scalar doesn't do anything, so we need to mark
1079 ! individual nested bcs to the du,dv,dw, whereas the vel_res can
1080 ! just get symmetry as a whole, because on this list we call
1081 ! apply_vector.
1082 ! Additionally we have to mark the special surface bc for p.
1083 call this%bclst_vel_res%append(bc_i)
1084 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1085 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1086 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1087
1088 call this%bcs_vel%append(bc_i)
1089
1090 call this%bc_sym_surface%mark_facets(bc_i%marked_facet)
1091 type is (non_normal_t)
1092 ! This is a bc for the residuals and increments, not the
1093 ! velocity itself. So, don't append to bcs_vel
1094 call this%bclst_vel_res%append(bc_i)
1095 call this%bc_du%mark_facets(bc_i%bc_x%marked_facet)
1096 call this%bc_dv%mark_facets(bc_i%bc_y%marked_facet)
1097 call this%bc_dw%mark_facets(bc_i%bc_z%marked_facet)
1098 type is (shear_stress_t)
1099 ! Same as symmetry
1100 call this%bclst_vel_res%append(bc_i%symmetry)
1101 call this%bclst_du%append(bc_i%symmetry%bc_x)
1102 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1103 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1104
1105 call this%bcs_vel%append(bc_i)
1106 type is (wall_model_bc_t)
1107 ! Same as symmetry
1108 call this%bclst_vel_res%append(bc_i%symmetry)
1109 call this%bclst_du%append(bc_i%symmetry%bc_x)
1110 call this%bclst_dv%append(bc_i%symmetry%bc_y)
1111 call this%bclst_dw%append(bc_i%symmetry%bc_z)
1112
1113 call this%bcs_vel%append(bc_i)
1114 class default
1115
1116 ! For the default case we use our dummy zero_dirichlet bcs to
1117 ! mark the same faces as in ordinary velocity dirichlet
1118 ! conditions.
1119 ! Additionally we mark the special PnPn pressure bc.
1120 if (bc_i%strong) then
1121 call this%bc_vel_res%mark_facets(bc_i%marked_facet)
1122 call this%bc_du%mark_facets(bc_i%marked_facet)
1123 call this%bc_dv%mark_facets(bc_i%marked_facet)
1124 call this%bc_dw%mark_facets(bc_i%marked_facet)
1125
1126 call this%bc_prs_surface%mark_facets(bc_i%marked_facet)
1127 end if
1128
1129 call this%bcs_vel%append(bc_i)
1130 end select
1131 end if
1132 end do
1133
1134 if (this%ale%active .and. (.not. this%ale%has_moving_boundary)) then
1135 call neko_error("Case file error: ALE is active, " // &
1136 "but no moving wall was found. " // &
1137 "Use type = 'no_slip' with 'moving': true in case file.")
1138 end if
1139
1140 ! Make sure all labeled zones with non-zero size have been marked
1141 do i = 1, size(this%msh%labeled_zones)
1142 if ((this%msh%labeled_zones(i)%size .gt. 0) .and. &
1143 (.not. marked_zones(i))) then
1144 write(error_unit, '(A, A, I0)') "*** ERROR ***: ", &
1145 "No fluid boundary condition assigned to zone ", i
1146 error stop
1147 end if
1148 end do
1149
1150 !
1151 ! Pressure bcs
1152 !
1153 call this%bcs_prs%init(n_bcs)
1154
1155 do i = 1, n_bcs
1156 ! Create a new json containing just the subdict for this bc
1157 call json_extract_item(core, bc_object, i, bc_subdict)
1158 bc_i => null()
1159 call pressure_bc_factory(bc_i, this, bc_subdict, this%c_Xh, user)
1160
1161 ! Not all bcs require an allocation for pressure in particular,
1162 ! so we check.
1163 if (associated(bc_i)) then
1164 call this%bcs_prs%append(bc_i)
1165
1166 ! Mark strong bcs in the dummy dp bc to force zero change.
1167 if (bc_i%strong) then
1168 call this%bc_dp%mark_facets(bc_i%marked_facet)
1169 end if
1170
1171 end if
1172
1173 end do
1174 else
1175 ! Check that there are no labeled zones, i.e. all are periodic.
1176 do i = 1, size(this%msh%labeled_zones)
1177 if (this%msh%labeled_zones(i)%size .gt. 0) then
1178 call neko_error("No boundary_conditions entry in the case file!")
1179 end if
1180 end do
1181
1182 ! For a pure periodic case, we still need to initilise the bc lists
1183 ! to a zero size to avoid issues with apply() in step()
1184 call this%bcs_vel%init()
1185 call this%bcs_prs%init()
1186
1187 end if
1188
1189 call this%bc_prs_surface%finalize()
1190 call this%bc_sym_surface%finalize()
1191
1192 call this%bc_vel_res%finalize()
1193 call this%bc_du%finalize()
1194 call this%bc_dv%finalize()
1195 call this%bc_dw%finalize()
1196 call this%bc_dp%finalize()
1197
1198 call this%bclst_vel_res%append(this%bc_vel_res)
1199 call this%bclst_du%append(this%bc_du)
1200 call this%bclst_dv%append(this%bc_dv)
1201 call this%bclst_dw%append(this%bc_dw)
1202 call this%bclst_dp%append(this%bc_dp)
1203
1204 ! If we have no strong pressure bcs, we will demean the pressure
1205 this%prs_dirichlet = .not. this%bclst_dp%is_empty()
1206 call mpi_allreduce(mpi_in_place, this%prs_dirichlet, 1, &
1207 mpi_logical, mpi_lor, neko_comm)
1208
1209
1210 if (allocated(marked_zones)) then
1211 deallocate(marked_zones)
1212 end if
1213
1214 if (allocated(zone_indices)) then
1215 deallocate(zone_indices)
1216 end if
1217
1218 end subroutine fluid_pnpn_setup_bcs
1219
1221 subroutine fluid_pnpn_write_boundary_conditions(this)
1222 use inflow, only : inflow_t
1224 use blasius, only : blasius_t
1226 use dong_outflow, only : dong_outflow_t
1227 use no_slip, only : no_slip_t
1228 class(fluid_pnpn_t), target, intent(inout) :: this
1229 type(dirichlet_t) :: bdry_mask
1230 type(field_t), pointer :: bdry_field
1231 type(file_t) :: bdry_file
1232 integer :: temp_index, i
1233 class(bc_t), pointer :: bci
1234 character(len=LOG_SIZE) :: log_buf
1235
1236 write(log_buf, '(A)') 'Marking using integer keys in bdry0.f00000'
1237 call neko_log%message(log_buf)
1238 write(log_buf, '(A)') 'Condition-value pairs: '
1239 call neko_log%message(log_buf)
1240 write(log_buf, '(A)') ' no_slip (stationary wall) = 1'
1241 call neko_log%message(log_buf)
1242 write(log_buf, '(A)') ' velocity_value = 2'
1243 call neko_log%message(log_buf)
1244 write(log_buf, '(A)') ' outflow, normal_outflow (+dong) = 3'
1245 call neko_log%message(log_buf)
1246 write(log_buf, '(A)') ' symmetry = 4'
1247 call neko_log%message(log_buf)
1248 write(log_buf, '(A)') ' periodic = 6'
1249 call neko_log%message(log_buf)
1250 write(log_buf, '(A)') ' user_velocity = 7'
1251 call neko_log%message(log_buf)
1252 write(log_buf, '(A)') ' user_pressure = 8'
1253 call neko_log%message(log_buf)
1254 write(log_buf, '(A)') ' shear_stress = 9'
1255 call neko_log%message(log_buf)
1256 write(log_buf, '(A)') ' wall_modelling = 10'
1257 call neko_log%message(log_buf)
1258 write(log_buf, '(A)') ' blasius_profile = 11'
1259 call neko_log%message(log_buf)
1260 write(log_buf, '(A)') ' no_slip (moving wall) = 12'
1261 call neko_log%message(log_buf)
1262
1263 call neko_scratch_registry%request_field(bdry_field, temp_index, .true.)
1264
1265
1266
1267 call bdry_mask%init_from_components(this%c_Xh, 6.0_rp)
1268 call bdry_mask%mark_zone(this%msh%periodic)
1269 call bdry_mask%finalize()
1270 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1271 call bdry_mask%free()
1272
1273 do i = 1, this%bcs_prs%size()
1274 bci => this%bcs_prs%get(i)
1275 select type (bc => bci)
1276 type is (zero_dirichlet_t)
1277 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1278 call bdry_mask%mark_facets(bci%marked_facet)
1279 call bdry_mask%finalize()
1280 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1281 call bdry_mask%free()
1282 type is (dong_outflow_t)
1283 call bdry_mask%init_from_components(this%c_Xh, 3.0_rp)
1284 call bdry_mask%mark_facets(bci%marked_facet)
1285 call bdry_mask%finalize()
1286 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1287 call bdry_mask%free()
1288 type is (field_dirichlet_t)
1289 call bdry_mask%init_from_components(this%c_Xh, 8.0_rp)
1290 call bdry_mask%mark_facets(bci%marked_facet)
1291 call bdry_mask%finalize()
1292 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1293 call bdry_mask%free()
1294 end select
1295 end do
1296
1297 do i = 1, this%bcs_vel%size()
1298 bci => this%bcs_vel%get(i)
1299 select type (bc => bci)
1300 type is (no_slip_t)
1301 if (bc%is_moving) then
1302 ! moving wall
1303 call bdry_mask%init_from_components(this%c_Xh, 12.0_rp)
1304 else
1305 ! stationary wall
1306 call bdry_mask%init_from_components(this%c_Xh, 1.0_rp)
1307 end if
1308 call bdry_mask%mark_facets(bci%marked_facet)
1309 call bdry_mask%finalize()
1310 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1311 call bdry_mask%free()
1312 type is (inflow_t)
1313 call bdry_mask%init_from_components(this%c_Xh, 2.0_rp)
1314 call bdry_mask%mark_facets(bci%marked_facet)
1315 call bdry_mask%finalize()
1316 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1317 call bdry_mask%free()
1318 type is (symmetry_t)
1319 call bdry_mask%init_from_components(this%c_Xh, 4.0_rp)
1320 call bdry_mask%mark_facets(bci%marked_facet)
1321 call bdry_mask%finalize()
1322 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1323 call bdry_mask%free()
1324 type is (field_dirichlet_vector_t)
1325 call bdry_mask%init_from_components(this%c_Xh, 7.0_rp)
1326 call bdry_mask%mark_facets(bci%marked_facet)
1327 call bdry_mask%finalize()
1328 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1329 call bdry_mask%free()
1330 type is (shear_stress_t)
1331 call bdry_mask%init_from_components(this%c_Xh, 9.0_rp)
1332 call bdry_mask%mark_facets(bci%marked_facet)
1333 call bdry_mask%finalize()
1334 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1335 call bdry_mask%free()
1336 type is (wall_model_bc_t)
1337 call bdry_mask%init_from_components(this%c_Xh, 10.0_rp)
1338 call bdry_mask%mark_facets(bci%marked_facet)
1339 call bdry_mask%finalize()
1340 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1341 call bdry_mask%free()
1342 type is (blasius_t)
1343 call bdry_mask%init_from_components(this%c_Xh, 11.0_rp)
1344 call bdry_mask%mark_facets(bci%marked_facet)
1345 call bdry_mask%finalize()
1346 call bdry_mask%apply_scalar(bdry_field%x, this%dm_Xh%size())
1347 call bdry_mask%free()
1348 end select
1349 end do
1350
1351
1352 call bdry_file%init('bdry.fld')
1353 call bdry_file%write(bdry_field)
1354
1355 call neko_scratch_registry%relinquish_field(temp_index)
1356 end subroutine fluid_pnpn_write_boundary_conditions
1357
1358end module fluid_pnpn
Copy data between host and device (or device and device)
Definition device.F90:72
Retrieves a parameter by name or assigns a provided default value. In the latter case also adds the m...
Retrieves a parameter by name or throws an error.
Apply cyclic boundary condition to a vector field.
Subroutines to add advection terms to the RHS of a transport equation.
Definition advection.f90:34
ALE Manager: Handles Mesh Motion.
Defines a Matrix-vector product.
Definition ax.f90:34
Defines a list of bc_t.
Definition bc_list.f90:34
Defines a boundary condition.
Definition bc.f90:34
Defines a Blasius profile dirichlet condition.
Definition blasius.f90:34
Defines a checkpoint.
Coefficients.
Definition coef.f90:34
Definition comm.F90:1
type(mpi_comm), public neko_comm
MPI communicator.
Definition comm.F90:45
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
subroutine, public device_event_sync(event)
Synchronize an event.
Definition device.F90:1594
integer, parameter, public host_to_device
Definition device.F90:48
type(c_ptr), bind(C), public glb_cmd_event
Event for the global command queue.
Definition device.F90:63
Defines a dirichlet boundary condition.
Definition dirichlet.f90:34
Defines a dong outflow condition.
Dirichlet condition applied in the facet normal direction.
Defines inflow dirichlet conditions.
Defines user dirichlet condition for a scalar field.
subroutine, public field_add2(a, b, n)
Vector addition .
subroutine, public field_copy(a, b, n)
Copy a vector .
Contains the field_serties_t type.
Defines a field.
Definition field.f90:34
Module for file I/O operations.
Definition file.f90:34
Auxiliary routines for fluid solvers.
Definition fluid_aux.f90:34
subroutine, public fluid_step_info(time, ksp_results, full_stress_formulation, strict_convergence, allow_stabilization, iteration)
Prints for prs, velx, vely, velz the following: Number of iterations, start residual,...
Definition fluid_aux.f90:54
Modular version of the Classic Nek5000 Pn/Pn formulation for fluids.
subroutine fluid_pnpn_setup_bcs(this, user, params)
Sets up the boundary condition for the scheme.
subroutine fluid_pnpn_restart(this, chkp)
subroutine fluid_pnpn_write_boundary_conditions(this)
Write a field with boundary condition specifications.
subroutine fluid_pnpn_step(this, time, dt_controller)
Advance fluid simulation in time.
subroutine fluid_pnpn_free(this)
subroutine fluid_pnpn_init(this, msh, lx, params, user, chkp)
Boundary condition factory for pressure.
Defines Gather-scatter operations.
Definition gs_ops.f90:34
integer, parameter, public gs_op_add
Definition gs_ops.f90:36
Defines inflow dirichlet conditions.
Definition inflow.f90:34
Utilities for retrieving parameters from the case files.
Implements the base abstract type for Krylov solvers plus helper types.
Definition krylov.f90:34
Logging routines.
Definition log.f90:34
type(log_t), public neko_log
Global log stream.
Definition log.f90:80
integer, parameter, public log_size
Definition log.f90:46
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:101
subroutine, public opadd2cm(a1, a2, a3, b1, b2, b3, c, n, gdim)
Definition mathops.f90:154
Defines a mesh.
Definition mesh.f90:34
Build configurations.
integer, parameter neko_bcknd_device
Defines no-slip boundary condition (extends zero_dirichlet)
Definition no_slip.f90:34
Dirichlet condition on axis aligned plane in the non normal direction.
integer, parameter, public dp
Definition num_types.f90:9
integer, parameter, public rp
Global precision used in computations.
Definition num_types.f90:12
Operators.
Definition operators.f90:34
subroutine, public ortho(x, glb_n_points, n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
Operators accelerator backends.
subroutine, public device_ortho(x_d, glb_n_points, n)
Othogonalize with regard to vector (1,1,1,1,1,1...,1)^T.
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:79
subroutine, public profiler_end_region(name, region_id)
End the most recently started profiler region.
Definition profiler.F90:116
Project x onto X , the space of old solutions and back again Couple projections for velocity.
Project x onto X, the space of old solutions and back again.
Defines a registry for storing solution fields.
Definition registry.f90:34
type(registry_t), target, public neko_registry
Global field registry.
Definition registry.f90:144
Routines to generate the right-hand sides for the convection-diffusion equation. Employs the EXT/BDF ...
Definition rhs_maker.f90:38
Defines a registry for storing and requesting temporary objects This can be used when you have a func...
type(scratch_registry_t), target, public neko_scratch_registry
Global scratch registry.
Defines a shear stress boundary condition for a vector field. Maintainer: Timofey Mukha.
Implements the source_term_t type and a wrapper source_term_wrapper_t.
Mixed Dirichlet-Neumann axis aligned symmetry plane.
Definition symmetry.f90:34
Module with things related to the simulation time.
Implements type time_step_controller.
Interfaces for user interaction with NEKO.
Definition user_intf.f90:34
Utilities.
Definition utils.f90:35
subroutine, public neko_type_error(base_type, wrong_type, known_types)
Reports an error allocating a type for a particular base pointer class.
Definition utils.f90:365
Defines the wall_model_bc_t type. Maintainer: Timofey Mukha.
Defines a zero-valued Dirichlet boundary condition.
Base abstract type for computing the advection operator.
Definition advection.f90:46
Base type for a matrix-vector product providing .
Definition ax.f90:43
Base type for a boundary condition.
Definition bc.f90:62
A list of allocatable `bc_t`. Follows the standard interface of lists.
Definition bc_list.f90:49
Blasius profile for inlet (vector valued).
Definition blasius.f90:55
Coefficients defined on a given (mesh, ) tuple. Arrays use indices (i,j,k,e): element e,...
Definition coef.f90:63
Generic Dirichlet boundary condition on .
Definition dirichlet.f90:49
Dong outflow condition Follows "A Convective-like Energy-Stable Open Boundary Condition for Simulati...
Dirichlet condition in facet normal direction.
User defined dirichlet condition, for which the user can work with an entire field....
Extension of the user defined dirichlet condition field_dirichlet
Stores a series (sequence) of fields, logically connected to a base field, and arranged according to ...
A wrapper around a polymorphic generic_file_t that handles its init. This is essentially a factory fo...
Definition file.f90:56
Dirichlet condition for inlet (vector valued)
Definition inflow.f90:47
Type for storing initial and final residuals in a Krylov solver.
Definition krylov.f90:56
Dirichlet condition in non normal direction of a plane.
Abstract type to compute pressure residual.
Definition pnpn_res.f90:48
Abstract type to compute velocity residual.
Definition pnpn_res.f90:54
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
A shear stress boundary condition.
Mixed Dirichlet-Neumann symmetry plane condition.
Definition symmetry.f90:50
A struct that contains all info about the time, expand as needed.
A type collecting all the overridable user routines and flag to suppress type injection from custom m...
A shear stress boundary condition, computing the stress values using a wall model.
Zero-valued Dirichlet boundary condition. Used for no-slip walls, but also for various auxillary cond...