Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_output.f90
Go to the documentation of this file.
1
34!
37 use num_types, only: rp
38 use adjoint_fluid_scheme, only: adjoint_fluid_scheme_t
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
47 implicit none
48 private
49
51 type, public, extends(output_t) :: adjoint_output_t
52 type(field_list_t) :: adjoint
53 logical :: always_write_mesh = .false.
54 contains
56 procedure, pass(this) :: init => adjoint_output_init
58 procedure, pass(this) :: free => adjoint_output_free
60 procedure, pass(this) :: sample => adjoint_output_sample
61 end type adjoint_output_t
62
63contains
64
76 subroutine adjoint_output_init(this, precision, adjoint, adjoint_scalars, &
77 name, path, fmt, layout, always_write_mesh)
78 class(adjoint_output_t), intent(inout) :: this
79 integer, intent(in) :: precision
80 class(adjoint_fluid_scheme_t), intent(in), target :: adjoint
81 class(adjoint_scalars_t), intent(in), optional, target :: adjoint_scalars
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
89
90 call this%free()
91
92 suffix = '.fld'
93 if (present(fmt)) then
94 if (fmt .eq. 'adios2') then
95 suffix = '.bp'
96 else if (fmt .eq. 'vtkhdf') then
97 suffix = '.vtkhdf'
98 end if
99 end if
100
101 if (present(always_write_mesh)) then
102 this%always_write_mesh = always_write_mesh
103 end if
104
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)
111 else
112 fname = 'field' // trim(suffix)
113 end if
114
115 if (present(layout)) then
116 call this%init_base(fname, precision, layout)
117 else
118 call this%init_base(fname, precision)
119 end if
120
121 ! Calculate total number of fields
122 n_scalars = 0
123 if (present(adjoint_scalars)) then
124 n_scalars = size(adjoint_scalars%adjoint_scalar_fields)
125 end if
126
127 ! Initialize field list with appropriate size
128 call this%adjoint%init(4 + n_scalars)
129
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)
134
135 ! Assign all scalar fields
136 if (present(adjoint_scalars)) then
137 do i = 1, n_scalars
138 call this%adjoint%assign(4 + i, &
139 adjoint_scalars%adjoint_scalar_fields(i)%s_adj)
140 end do
141 end if
142
143 end subroutine adjoint_output_init
144
146 subroutine adjoint_output_free(this)
147 class(adjoint_output_t), intent(inout) :: this
148
149 call this%free_base()
150 call this%adjoint%free()
151
152 end subroutine adjoint_output_free
153
157 subroutine adjoint_output_sample(this, t)
158 class(adjoint_output_t), intent(inout) :: this
159 real(kind=rp), intent(in) :: t
160 integer :: i
161
162 if (neko_bcknd_device .eq. 1) then
163
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))) ! Sync on the last field
169 end do
170 end associate
171
172 end if
173
174 select type (ft => this%file_%file_type)
175 ! Only fld files have the option to write the mesh at command
176 type is (fld_file_t)
177 ft%write_mesh = this%always_write_mesh
178 call ft%write(this%adjoint, t)
179 class default
180 call ft%write(this%adjoint, t)
181 end select
182
183 end subroutine adjoint_output_sample
184
185end module adjoint_output
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 for a scalar advection-diffusion solver.
Type to manage multiple adjoint scalar transport equations.