36  use gather_scatter, 
only: gs_t, gs_op_min, gs_op_max
 
   37  use neko_config, 
only: neko_bcknd_device
 
   38  use num_types, 
only: rp, i8
 
   40  use field, 
only: field_t
 
   41  use space, 
only: space_t, gll
 
   42  use dofmap, 
only: dofmap_t
 
   43  use krylov, 
only: ksp_t, krylov_solver_factory, ksp_max_iter
 
   44  use coefs, 
only: coef_t
 
   45  use jacobi, 
only: jacobi_t
 
   46  use sx_jacobi, 
only: sx_jacobi_t
 
   47  use device_jacobi, 
only: device_jacobi_t
 
   48  use hsmg, 
only: hsmg_t
 
   49  use phmg, 
only: phmg_t
 
   50  use precon, 
only: pc_t, precon_factory, precon_destroy
 
   51  use fluid_stats, 
only: fluid_stats_t
 
   53  use bc_list, 
only: bc_list_t
 
   54  use mesh, 
only: mesh_t, neko_msh_max_zlbls, neko_msh_max_zlbl_len
 
   55  use operators, 
only: cfl
 
   56  use logger, 
only: neko_log, log_size, neko_log_verbose
 
   57  use field_registry, 
only: neko_field_registry
 
   58  use json_utils, 
only: json_get, json_get_or_default
 
   59  use json_module, 
only: json_file
 
   60  use scratch_registry, 
only: scratch_registry_t
 
   61  use user_intf, 
only: user_t, dummy_user_material_properties, &
 
   62       user_material_properties_intf
 
   63  use utils, 
only: neko_error
 
   64  use time_state, 
only: time_state_t
 
   67  use device_math, 
only: device_cfill, device_add2s2
 
   68  use field_math, 
only: field_cfill, field_add2s2, field_addcol3
 
   71  use device, 
only : device_event_sync, glb_cmd_event, device_to_host, &
 
   81     class(ksp_t), 
allocatable :: ksp_vel
 
   82     class(ksp_t), 
allocatable :: ksp_prs
 
   83     class(pc_t), 
allocatable :: pc_vel
 
   84     class(pc_t), 
allocatable :: pc_prs
 
   85     integer :: vel_projection_dim
 
   86     integer :: pr_projection_dim
 
   87     integer :: vel_projection_activ_step
 
   88     integer :: pr_projection_activ_step
 
   89     logical :: strict_convergence
 
   92     type(field_t), 
pointer :: u_adj_e => null() 
 
   93     type(field_t), 
pointer :: v_adj_e => null() 
 
   94     type(field_t), 
pointer :: w_adj_e => null() 
 
   96     type(fluid_stats_t) :: stats
 
   97     logical :: forced_flow_rate = .false. 
 
  100     character(len=:), 
allocatable :: nut_field_name
 
  103     integer(kind=i8) :: glb_n_points
 
  105     integer(kind=i8) :: glb_unique_points
 
  106     type(scratch_registry_t) :: scratch
 
  109     procedure, pass(this) :: init_base => adjoint_fluid_scheme_init_base
 
  110     procedure, pass(this) :: scheme_free => adjoint_fluid_scheme_free
 
  112     procedure, pass(this) :: validate => adjoint_fluid_scheme_validate
 
  114     procedure, pass(this) :: bc_apply_vel => adjoint_fluid_scheme_bc_apply_vel
 
  116     procedure, pass(this) :: bc_apply_prs => adjoint_fluid_scheme_bc_apply_prs
 
  118     procedure, pass(this) :: compute_cfl => adjoint_compute_cfl
 
  120     procedure, pass(this) :: set_material_properties => &
 
  121          adjoint_fluid_scheme_set_material_properties
 
  124     procedure, pass(this) :: update_material_properties => &
 
  125          adjoint_fluid_scheme_update_material_properties
 
  127     procedure, 
nopass :: solver_factory => adjoint_fluid_scheme_solver_factory
 
  129     procedure, pass(this) :: precon_factory_ => &
 
  130          adjoint_fluid_scheme_precon_factory
 
 
  135     module subroutine adjoint_fluid_scheme_factory(object, type_name)
 
  137       character(len=*) :: type_name
 
  138     end subroutine adjoint_fluid_scheme_factory
 
  141  public :: adjoint_fluid_scheme_incompressible_t, adjoint_fluid_scheme_factory
 
  146  subroutine adjoint_fluid_scheme_init_base(this, msh, lx, params, scheme, &
 
  148    class(adjoint_fluid_scheme_incompressible_t), 
target, 
intent(inout) :: this
 
  149    type(mesh_t), 
target, 
intent(inout) :: msh
 
  150    integer, 
intent(in) :: lx
 
  151    character(len=*), 
intent(in) :: scheme
 
  152    type(json_file), 
target, 
intent(inout) :: params
 
  153    type(user_t), 
target, 
intent(in) :: user
 
  154    logical, 
intent(in) :: kspv_init
 
  155    character(len=LOG_SIZE) :: log_buf
 
  156    real(kind=rp) :: real_val
 
  157    logical :: logical_val
 
  158    integer :: integer_val
 
  159    character(len=:), 
allocatable :: string_val1, string_val2
 
  160    character(len=:), 
allocatable :: json_key
 
  161    type(json_file) :: json_subdict
 
  169    if (msh%gdim .eq. 2) 
then 
  170       call this%Xh%init(gll, lx, lx)
 
  172       call this%Xh%init(gll, lx, lx, lx)
 
  175    call this%dm_Xh%init(msh, this%Xh)
 
  177    call this%gs_Xh%init(this%dm_Xh)
 
  179    call this%c_Xh%init(this%gs_Xh)
 
  182    call this%scratch%init(this%dm_Xh, 10, 2)
 
  185    call json_get_or_default(params, 
'case.fluid.name', this%name, 
"fluid")
 
  191    call neko_log%section(
'Adjoint fluid')
 
  192    write(log_buf, 
'(A, A)') 
'Type       : ', trim(scheme)
 
  193    call neko_log%message(log_buf)
 
  194    write(log_buf, 
'(A, A)') 
'Name       : ', trim(this%name)
 
  195    call neko_log%message(log_buf)
 
  198    call neko_field_registry%add_field(this%dm_Xh, 
'u_adj')
 
  199    call neko_field_registry%add_field(this%dm_Xh, 
'v_adj')
 
  200    call neko_field_registry%add_field(this%dm_Xh, 
'w_adj')
 
  201    this%u_adj => neko_field_registry%get_field(
'u_adj')
 
  202    this%v_adj => neko_field_registry%get_field(
'v_adj')
 
  203    this%w_adj => neko_field_registry%get_field(
'w_adj')
 
  208    call this%set_material_properties(params, user)
 
  212         'case.adjoint_fluid.velocity_solver.projection_space_size', &
 
  213         'case.fluid.velocity_solver.projection_space_size')
 
  214    call json_get_or_default(params, json_key, this%vel_projection_dim, 0)
 
  216         'case.adjoint_fluid.pressure_solver.projection_space_size', &
 
  217         'case.fluid.pressure_solver.projection_space_size')
 
  218    call json_get_or_default(params, json_key, this%pr_projection_dim, 0)
 
  221         'case.adjoint_fluid.velocity_solver.projection_hold_steps', &
 
  222         'case.fluid.velocity_solver.projection_hold_steps')
 
  223    call json_get_or_default(params, json_key, &
 
  224         this%vel_projection_activ_step, 5)
 
  226         'case.adjoint_fluid.pressure_solver.projection_hold_steps', &
 
  227         'case.fluid.pressure_solver.projection_hold_steps')
 
  228    call json_get_or_default(params, json_key, &
 
  229         this%pr_projection_activ_step, 5)
 
  233    call json_get_or_default(params, json_key, this%freeze, .false.)
 
  242    if (params%valid_path(
"case.fluid.flow_rate_force")) 
then 
  243       call neko_error(
'Flow rate forcing not yet implemented')
 
  244       this%forced_flow_rate = .true.
 
  249       write(log_buf, 
'(A, I1)') 
'Poly order : ', lx-1
 
  250    else if (lx .ge. 10) 
then 
  251       write(log_buf, 
'(A, I2)') 
'Poly order : ', lx-1
 
  253       write(log_buf, 
'(A, I3)') 
'Poly order : ', lx-1
 
  255    call neko_log%message(log_buf)
 
  256    this%glb_n_points = int(this%msh%glb_nelv, i8)*int(this%Xh%lxyz, i8)
 
  257    this%glb_unique_points = int(glsum(this%c_Xh%mult, this%dm_Xh%size()), i8)
 
  259    write(log_buf, 
'(A, I0)') 
'GLL points : ', this%glb_n_points
 
  260    call neko_log%message(log_buf)
 
  261    write(log_buf, 
'(A, I0)') 
'Unique pts.: ', this%glb_unique_points
 
  262    call neko_log%message(log_buf)
 
  264    call json_get(params, 
'case.numerics.dealias', logical_val)
 
  265    write(log_buf, 
'(A, L1)') 
'Dealias    : ', logical_val
 
  266    call neko_log%message(log_buf)
 
  268    call json_get_or_default(params, 
'case.output_boundary', logical_val, &
 
  270    write(log_buf, 
'(A, L1)') 
'Save bdry  : ', logical_val
 
  271    call neko_log%message(log_buf)
 
  273    call json_get_or_default(params, 
"case.fluid.full_stress_formulation", &
 
  274         logical_val, .false.)
 
  275    write(log_buf, 
'(A, L1)') 
'Full stress: ', logical_val
 
  276    call neko_log%message(log_buf)
 
  281    allocate(this%f_adj_x)
 
  282    allocate(this%f_adj_y)
 
  283    allocate(this%f_adj_z)
 
  284    call this%f_adj_x%init(this%dm_Xh, fld_name = 
"adjoint_rhs_x")
 
  285    call this%f_adj_y%init(this%dm_Xh, fld_name = 
"adjoint_rhs_y")
 
  286    call this%f_adj_z%init(this%dm_Xh, fld_name = 
"adjoint_rhs_z")
 
  306       call neko_log%section(
"Adjoint Velocity solver")
 
  309            'case.adjoint_fluid.velocity_solver.max_iterations', &
 
  310            'case.fluid.velocity_solver.max_iterations')
 
  311       call json_get_or_default(params, json_key, integer_val, ksp_max_iter)
 
  314            'case.adjoint_fluid.velocity_solver.type', &
 
  315            'case.fluid.velocity_solver.type')
 
  316       call json_get(params, json_key, string_val1)
 
  319            'case.adjoint_fluid.velocity_solver.preconditioner.type', &
 
  320            'case.fluid.velocity_solver.preconditioner.type')
 
  321       call json_get(params, json_key, string_val2)
 
  324            'case.adjoint_fluid.velocity_solver.preconditioner', &
 
  325            'case.fluid.velocity_solver.preconditioner')
 
  326       call json_get(params, json_key, json_subdict)
 
  329            'case.adjoint_fluid.velocity_solver.absolute_tolerance', &
 
  330            'case.fluid.velocity_solver.absolute_tolerance')
 
  331       call json_get(params, json_key, real_val)
 
  334            'case.adjoint_fluid.velocity_solver.monitor', &
 
  335            'case.fluid.velocity_solver.monitor')
 
  336       call json_get_or_default(params, json_key, logical_val, .false.)
 
  338       call neko_log%message(
'Type       : (' // trim(string_val1) // &
 
  339            ', ' // trim(string_val2) // 
')')
 
  341       write(log_buf, 
'(A,ES13.6)') 
'Abs tol    :', real_val
 
  342       call neko_log%message(log_buf)
 
  343       call this%solver_factory(this%ksp_vel, this%dm_Xh%size(), &
 
  344            string_val1, integer_val, real_val, logical_val)
 
  345       call this%precon_factory_(this%pc_vel, this%ksp_vel, &
 
  346            this%c_Xh, this%dm_Xh, this%gs_Xh, this%bcs_vel, &
 
  347            string_val2, json_subdict)
 
  348       call neko_log%end_section()
 
  352    call json_get_or_default(params, 
'case.fluid.strict_convergence', &
 
  353         this%strict_convergence, .false.)
 
  356    call this%ulag%init(this%u_adj, 2)
 
  357    call this%vlag%init(this%v_adj, 2)
 
  358    call this%wlag%init(this%w_adj, 2)
 
  360    call neko_field_registry%add_field(this%dm_Xh, 
'u_adj_e')
 
  361    call neko_field_registry%add_field(this%dm_Xh, 
'v_adj_e')
 
  362    call neko_field_registry%add_field(this%dm_Xh, 
'w_adj_e')
 
  363    this%u_adj_e => neko_field_registry%get_field(
'u_adj_e')
 
  364    this%v_adj_e => neko_field_registry%get_field(
'v_adj_e')
 
  365    this%w_adj_e => neko_field_registry%get_field(
'w_adj_e')
 
  368    call this%source_term%init(this%f_adj_x, this%f_adj_y, this%f_adj_z, &
 
  369         this%c_Xh, user, this%name)
 
  370    call this%source_term%add(params, 
'case.adjoint_fluid.source_term')
 
  372    call neko_log%end_section()
 
 
  374  end subroutine adjoint_fluid_scheme_init_base
 
  377  subroutine adjoint_fluid_scheme_free(this)
 
  378    class(adjoint_fluid_scheme_incompressible_t), 
intent(inout) :: this
 
  382    if (
allocated(this%ksp_vel)) 
then 
  383       call this%ksp_vel%free()
 
  384       deallocate(this%ksp_vel)
 
  387    if (
allocated(this%ksp_prs)) 
then 
  388       call this%ksp_prs%free()
 
  389       deallocate(this%ksp_prs)
 
  392    if (
allocated(this%pc_vel)) 
then 
  393       call precon_destroy(this%pc_vel)
 
  394       deallocate(this%pc_vel)
 
  397    if (
allocated(this%pc_prs)) 
then 
  398       call precon_destroy(this%pc_prs)
 
  399       deallocate(this%pc_prs)
 
  402    call this%source_term%free()
 
  404    call this%gs_Xh%free()
 
  406    call this%c_Xh%free()
 
  408    call this%scratch%free()
 
  415    nullify(this%u_adj_e)
 
  416    nullify(this%v_adj_e)
 
  417    nullify(this%w_adj_e)
 
  419    call this%ulag%free()
 
  420    call this%vlag%free()
 
  421    call this%wlag%free()
 
  424    if (
associated(this%f_adj_x)) 
then 
  425       call this%f_adj_x%free()
 
  428    if (
associated(this%f_adj_y)) 
then 
  429       call this%f_adj_y%free()
 
  432    if (
associated(this%f_adj_z)) 
then 
  433       call this%f_adj_z%free()
 
  436    nullify(this%f_adj_x)
 
  437    nullify(this%f_adj_y)
 
  438    nullify(this%f_adj_z)
 
  443  end subroutine adjoint_fluid_scheme_free
 
  447  subroutine adjoint_fluid_scheme_validate(this)
 
  448    class(adjoint_fluid_scheme_incompressible_t), 
target, 
intent(inout) :: this
 
  451    if ((.not. 
associated(this%u_adj)) .or. &
 
  452         (.not. 
associated(this%v_adj)) .or. &
 
  453         (.not. 
associated(this%w_adj)) .or. &
 
  454         (.not. 
associated(this%p_adj))) 
then 
  455       call neko_error(
'Fields are not registered')
 
  458    if ((.not. 
allocated(this%u_adj%x)) .or. &
 
  459         (.not. 
allocated(this%v_adj%x)) .or. &
 
  460         (.not. 
allocated(this%w_adj%x)) .or. &
 
  461         (.not. 
allocated(this%p_adj%x))) 
then 
  462       call neko_error(
'Fields are not allocated')
 
  465    if (.not. 
allocated(this%ksp_vel)) 
then 
  466       call neko_error(
'No Krylov solver for velocity defined')
 
  469    if (.not. 
allocated(this%ksp_prs)) 
then 
  470       call neko_error(
'No Krylov solver for pressure defined')
 
  476    call this%chkp%init(this%u_adj, this%v_adj, this%w_adj, this%p_adj)
 
  478  end subroutine adjoint_fluid_scheme_validate
 
  484  subroutine adjoint_fluid_scheme_bc_apply_vel(this, time, strong)
 
  485    class(adjoint_fluid_scheme_incompressible_t), 
intent(inout) :: this
 
  486    type(time_state_t), 
intent(in) :: time
 
  487    logical, 
intent(in) :: strong
 
  490    class(bc_t), 
pointer :: b => null()
 
  492    call this%bcs_vel%apply_vector(&
 
  493         this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
 
  495    call this%gs_Xh%op(this%u_adj, gs_op_min, glb_cmd_event)
 
  496    call this%gs_Xh%op(this%v_adj, gs_op_min, glb_cmd_event)
 
  497    call this%gs_Xh%op(this%w_adj, gs_op_min, glb_cmd_event)
 
  498    call device_event_sync(glb_cmd_event)
 
  500    call this%bcs_vel%apply_vector(&
 
  501         this%u_adj%x, this%v_adj%x, this%w_adj%x, this%dm_Xh%size(), &
 
  503    call this%gs_Xh%op(this%u_adj, gs_op_max, glb_cmd_event)
 
  504    call this%gs_Xh%op(this%v_adj, gs_op_max, glb_cmd_event)
 
  505    call this%gs_Xh%op(this%w_adj, gs_op_max, glb_cmd_event)
 
  506    call device_event_sync(glb_cmd_event)
 
  508    do i = 1, this%bcs_vel%size()
 
  509       b => this%bcs_vel%get(i)
 
  515  end subroutine adjoint_fluid_scheme_bc_apply_vel
 
  519  subroutine adjoint_fluid_scheme_bc_apply_prs(this, time)
 
  520    class(adjoint_fluid_scheme_incompressible_t), 
intent(inout) :: this
 
  521    type(time_state_t), 
intent(in) :: time
 
  524    class(bc_t), 
pointer :: b => null()
 
  526    call this%bcs_prs%apply(this%p_adj, time)
 
  527    call this%gs_Xh%op(this%p_adj, gs_op_min, glb_cmd_event)
 
  528    call device_event_sync(glb_cmd_event)
 
  530    call this%bcs_prs%apply(this%p_adj, time)
 
  531    call this%gs_Xh%op(this%p_adj, gs_op_max, glb_cmd_event)
 
  532    call device_event_sync(glb_cmd_event)
 
  534    do i = 1, this%bcs_prs%size()
 
  535       b => this%bcs_prs%get(i)
 
  540  end subroutine adjoint_fluid_scheme_bc_apply_prs
 
  544  subroutine adjoint_fluid_scheme_solver_factory(ksp, n, solver, &
 
  545       max_iter, abstol, monitor)
 
  546    class(ksp_t), 
allocatable, 
target, 
intent(inout) :: ksp
 
  547    integer, 
intent(in), 
value :: n
 
  548    character(len=*), 
intent(in) :: solver
 
  549    integer, 
intent(in) :: max_iter
 
  550    real(kind=rp), 
intent(in) :: abstol
 
  551    logical, 
intent(in) :: monitor
 
  553    call krylov_solver_factory(ksp, n, solver, max_iter, abstol, &
 
  556  end subroutine adjoint_fluid_scheme_solver_factory
 
  559  subroutine adjoint_fluid_scheme_precon_factory(this, pc, ksp, coef, dof, gs, &
 
  560       bclst, pctype, pcparams)
 
  561    class(adjoint_fluid_scheme_incompressible_t), 
intent(inout) :: this
 
  562    class(pc_t), 
allocatable, 
target, 
intent(inout) :: pc
 
  563    class(ksp_t), 
target, 
intent(inout) :: ksp
 
  564    type(coef_t), 
target, 
intent(in) :: coef
 
  565    type(dofmap_t), 
target, 
intent(in) :: dof
 
  566    type(gs_t), 
target, 
intent(inout) :: gs
 
  567    type(bc_list_t), 
target, 
intent(inout) :: bclst
 
  568    character(len=*) :: pctype
 
  569    type(json_file), 
intent(inout) :: pcparams
 
  571    call precon_factory(pc, pctype)
 
  573    select type (pcp => pc)
 
  575       call pcp%init(coef, dof, gs)
 
  576    type is (sx_jacobi_t)
 
  577       call pcp%init(coef, dof, gs)
 
  578    type is (device_jacobi_t)
 
  579       call pcp%init(coef, dof, gs)
 
  581       call pcp%init(coef, bclst, pcparams)
 
  583       call pcp%init(coef, bclst, pcparams)
 
  588  end subroutine adjoint_fluid_scheme_precon_factory
 
  607  function adjoint_compute_cfl(this, dt) 
result(c)
 
  608    class(adjoint_fluid_scheme_incompressible_t), 
intent(in) :: this
 
  609    real(kind=rp), 
intent(in) :: dt
 
  612    c = cfl(dt, this%u_adj%x, this%v_adj%x, this%w_adj%x, &
 
  613         this%Xh, this%c_Xh, this%msh%nelv, this%msh%gdim)
 
  615  end function adjoint_compute_cfl
 
  621  subroutine adjoint_fluid_scheme_update_material_properties(this, time)
 
  622    class(adjoint_fluid_scheme_incompressible_t), 
intent(inout) :: this
 
  623    type(time_state_t), 
intent(in) :: time
 
  624    type(field_t), 
pointer :: nut
 
  626    call this%user_material_properties(this%name, &
 
  627         this%material_properties, time)
 
  629    if (len(trim(this%nut_field_name)) > 0) 
then 
  630       nut => neko_field_registry%get_field(this%nut_field_name)
 
  631       call field_addcol3(this%mu, nut, this%rho)
 
  638    if (neko_bcknd_device .eq. 1) 
then 
  639       call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
 
  640            device_to_host, sync=.false.)
 
  641       call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
 
  642            device_to_host, sync=.false.)
 
  644  end subroutine adjoint_fluid_scheme_update_material_properties
 
  650  subroutine adjoint_fluid_scheme_set_material_properties(this, params, user)
 
  651    class(adjoint_fluid_scheme_incompressible_t), 
intent(inout) :: this
 
  652    type(json_file), 
intent(inout) :: params
 
  653    type(user_t), 
target, 
intent(in) :: user
 
  654    character(len=LOG_SIZE) :: log_buf
 
  656    procedure(user_material_properties_intf), 
pointer :: dummy_mp_ptr
 
  657    real(kind=rp) :: const_mu, const_rho
 
  658    type(time_state_t) :: time
 
  661    dummy_mp_ptr => dummy_user_material_properties
 
  663    call this%mu%init(this%dm_Xh, 
"mu")
 
  664    call this%rho%init(this%dm_Xh, 
"rho")
 
  665    call this%material_properties%init(2)
 
  666    call this%material_properties%assign_to_field(1, this%rho)
 
  667    call this%material_properties%assign_to_field(2, this%mu)
 
  669    if (.not. 
associated(user%material_properties, dummy_mp_ptr)) 
then 
  671       write(log_buf, 
'(A)') 
"Material properties must be set in the user& 
  673       call neko_log%message(log_buf)
 
  674       this%user_material_properties => user%material_properties
 
  676       call user%material_properties(this%name, &
 
  677            this%material_properties, time)
 
  680       this%user_material_properties => dummy_user_material_properties
 
  682       if (params%valid_path(
'case.fluid.Re') .and. &
 
  683            (params%valid_path(
'case.fluid.mu') .or. &
 
  684            params%valid_path(
'case.fluid.rho'))) 
then 
  685          call neko_error(
"To set the material properties for the fluid, " // &
 
  686               "either provide Re OR mu and rho in the case file.")
 
  688       else if (params%valid_path(
'case.fluid.Re')) 
then 
  690          write(log_buf, 
'(A)') 
'Non-dimensional fluid material properties & 
  692          call neko_log%message(log_buf, lvl = neko_log_verbose)
 
  693          write(log_buf, 
'(A)') 
'Density will be set to 1, dynamic viscosity to& 
  695          call neko_log%message(log_buf, lvl = neko_log_verbose)
 
  698          call json_get(params, 
'case.fluid.Re', const_mu)
 
  699          write(log_buf, 
'(A)') 
'Read non-dimensional material properties' 
  700          call neko_log%message(log_buf)
 
  701          write(log_buf, 
'(A,ES13.6)') 
'Re         :', const_mu
 
  702          call neko_log%message(log_buf)
 
  707          const_mu = 1.0_rp/const_mu
 
  710          call json_get(params, 
'case.fluid.mu', const_mu)
 
  711          call json_get(params, 
'case.fluid.rho', const_rho)
 
  717    if (
associated(user%material_properties, dummy_mp_ptr)) 
then 
  719       call field_cfill(this%mu, const_mu)
 
  720       call field_cfill(this%rho, const_rho)
 
  723       write(log_buf, 
'(A,ES13.6)') 
'rho        :', const_rho
 
  724       call neko_log%message(log_buf)
 
  725       write(log_buf, 
'(A,ES13.6)') 
'mu         :', const_mu
 
  726       call neko_log%message(log_buf)
 
  733    if (neko_bcknd_device .eq. 1) 
then 
  734       call device_memcpy(this%rho%x, this%rho%x_d, this%rho%size(), &
 
  735            device_to_host, sync=.false.)
 
  736       call device_memcpy(this%mu%x, this%mu%x_d, this%mu%size(), &
 
  737            device_to_host, sync=.false.)
 
  739  end subroutine adjoint_fluid_scheme_set_material_properties
 
Implements the adjoint_source_term_t type.
 
Base type of all fluid formulations.
 
Base type of all fluid formulations.
 
Wrapper contaning and executing the adjoint source terms.