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