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 implicit none
47 private
48
50 type, public, extends(output_t) :: adjoint_output_t
51 type(field_list_t) :: adjoint
52 contains
54 procedure, pass(this) :: init => adjoint_output_init
56 procedure, pass(this) :: free => adjoint_output_free
58 procedure, pass(this) :: sample => adjoint_output_sample
59 end type adjoint_output_t
60
61contains
62
71 subroutine adjoint_output_init(this, precision, adjoint, adjoint_scalars, &
72 name, path)
73 class(adjoint_output_t), intent(inout) :: this
74 integer, intent(in) :: precision
75 class(adjoint_fluid_scheme_t), intent(in), target :: adjoint
76 class(adjoint_scalars_t), intent(in), optional, target :: adjoint_scalars
77 character(len=*), intent(in), optional :: name
78 character(len=*), intent(in), optional :: path
79 character(len=1024) :: fname
80 integer :: i, n_scalars
81
82 call this%free()
83
84 if (present(name) .and. present(path)) then
85 fname = trim(path) // trim(name) // '.fld'
86 else if (present(name)) then
87 fname = trim(name) // '.fld'
88 else if (present(path)) then
89 fname = trim(path) // 'adjoint.fld'
90 else
91 fname = 'adjoint.fld'
92 end if
93
94 call this%init_base(fname, precision)
95
96 ! Calculate total number of fields
97 n_scalars = 0
98 if (present(adjoint_scalars)) then
99 n_scalars = size(adjoint_scalars%adjoint_scalar_fields)
100 end if
101
102 ! Initialize field list with appropriate size
103 call this%adjoint%init(4 + n_scalars)
104
105 call this%adjoint%assign(1, adjoint%p_adj)
106 call this%adjoint%assign(2, adjoint%u_adj)
107 call this%adjoint%assign(3, adjoint%v_adj)
108 call this%adjoint%assign(4, adjoint%w_adj)
109
110 ! Assign all scalar fields
111 if (present(adjoint_scalars)) then
112 do i = 1, n_scalars
113 call this%adjoint%assign(4 + i, &
114 adjoint_scalars%adjoint_scalar_fields(i)%s_adj)
115 end do
116 end if
117
118 end subroutine adjoint_output_init
119
121 subroutine adjoint_output_free(this)
122 class(adjoint_output_t), intent(inout) :: this
123
124 call this%free_base()
125 call this%adjoint%free()
126
127 end subroutine adjoint_output_free
128
132 subroutine adjoint_output_sample(this, t)
133 class(adjoint_output_t), intent(inout) :: this
134 real(kind=rp), intent(in) :: t
135 integer :: i
136
137 if (neko_bcknd_device .eq. 1) then
138
139 associate(fields => this%adjoint%items)
140 do i = 1, size(fields)
141 call device_memcpy(fields(i)%ptr%x, fields(i)%ptr%x_d, &
142 fields(i)%ptr%dof%size(), device_to_host, &
143 sync = (i .eq. size(fields))) ! Sync on the last field
144 end do
145 end associate
146
147 end if
148
149 call this%file_%write(this%adjoint, t)
150
151 end subroutine adjoint_output_sample
152
153end module adjoint_output
Defines an output for a adjoint.
subroutine adjoint_output_init(this, precision, adjoint, adjoint_scalars, name, path)
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.