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