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