37 use num_types,
only: rp
39 use scalar_scheme,
only: scalar_scheme_t
41 use field_list,
only: field_list_t
42 use neko_config,
only: neko_bcknd_device
43 use device,
only: device_memcpy, device_to_host
44 use output,
only: output_t
46 use fld_file,
only: fld_file_t
52 type(field_list_t) :: adjoint
53 logical :: always_write_mesh = .false.
58 procedure, pass(this) :: free => adjoint_output_free
60 procedure, pass(this) :: sample => adjoint_output_sample
77 name, path, fmt, layout, always_write_mesh)
79 integer,
intent(in) :: precision
82 character(len=*),
intent(in),
optional :: name
83 character(len=*),
intent(in),
optional :: path
84 character(len=*),
intent(in),
optional :: fmt
85 logical,
intent(in),
optional :: always_write_mesh
86 integer,
intent(in),
optional :: layout
87 character(len=1024) :: fname, suffix
88 integer :: i, n_scalars
93 if (
present(fmt))
then
94 if (fmt .eq.
'adios2')
then
96 else if (fmt .eq.
'vtkhdf')
then
101 if (
present(always_write_mesh))
then
102 this%always_write_mesh = always_write_mesh
105 if (
present(name) .and.
present(path))
then
106 fname = trim(path) // trim(name) // trim(suffix)
107 else if (
present(name))
then
108 fname = trim(name) // trim(suffix)
109 else if (
present(path))
then
110 fname = trim(path) //
'field' // trim(suffix)
112 fname =
'field' // trim(suffix)
115 if (
present(layout))
then
116 call this%init_base(fname, precision, layout)
118 call this%init_base(fname, precision)
128 call this%adjoint%init(4 + n_scalars)
130 call this%adjoint%assign(1, adjoint%p_adj)
131 call this%adjoint%assign(2, adjoint%u_adj)
132 call this%adjoint%assign(3, adjoint%v_adj)
133 call this%adjoint%assign(4, adjoint%w_adj)
138 call this%adjoint%assign(4 + i, &
146 subroutine adjoint_output_free(this)
149 call this%free_base()
150 call this%adjoint%free()
152 end subroutine adjoint_output_free
157 subroutine adjoint_output_sample(this, t)
159 real(kind=rp),
intent(in) :: t
162 if (neko_bcknd_device .eq. 1)
then
164 associate(fields => this%adjoint%items)
165 do i = 1,
size(fields)
166 call device_memcpy(fields(i)%ptr%x, fields(i)%ptr%x_d, &
167 fields(i)%ptr%dof%size(), device_to_host, &
168 sync = (i .eq.
size(fields)))
174 select type (ft => this%file_%file_type)
177 ft%write_mesh = this%always_write_mesh
178 call ft%write(this%adjoint, t)
180 call ft%write(this%adjoint, t)
183 end subroutine adjoint_output_sample
Defines an output for a adjoint.
subroutine adjoint_output_init(this, precision, adjoint, adjoint_scalars, name, path, fmt, layout, always_write_mesh)
Constructor.
Contains the adjoint_scalar_scheme_t type.
Contains the adjoint_scalars_t type that manages multiple scalar fields.
Base type of all fluid formulations.
Base type for a scalar advection-diffusion solver.
Type to manage multiple adjoint scalar transport equations.