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