Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
adjoint_output.f90
1! Copyright (c) 2020-2023, The Neko Authors
2! All rights reserved.
3!
4! Redistribution and use in source and binary forms, with or without
5! modification, are permitted provided that the following conditions
6! are met:
7!
8! * Redistributions of source code must retain the above copyright
9! notice, this list of conditions and the following disclaimer.
10!
11! * Redistributions in binary form must reproduce the above
12! copyright notice, this list of conditions and the following
13! disclaimer in the documentation and/or other materials provided
14! with the distribution.
15!
16! * Neither the name of the authors nor the names of its
17! contributors may be used to endorse or promote products derived
18! from this software without specific prior written permission.
19!
20! THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS
21! "AS IS" AND ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT
22! LIMITED TO, THE IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS
23! FOR A PARTICULAR PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL THE
24! COPYRIGHT OWNER OR CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT,
25! INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING,
26! BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
27! LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER
28! CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
29! LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN
30! ANY WAY OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE
31! POSSIBILITY OF SUCH DAMAGE.
32!
35 use num_types, only: rp
36 use adjoint_fluid_scheme, only: adjoint_fluid_scheme_t
37 use scalar_scheme, only: scalar_scheme_t
39 use field_list, only: field_list_t
40 use neko_config, only: neko_bcknd_device
41 use device, only: device_memcpy, device_to_host
42 use output, only: output_t
44 implicit none
45 private
46
48 type, public, extends(output_t) :: adjoint_output_t
49 type(field_list_t) :: adjoint
50 contains
51 procedure, pass(this) :: sample => adjoint_output_sample
52 end type adjoint_output_t
53
54 interface adjoint_output_t
55 module procedure adjoint_output_init
56 end interface adjoint_output_t
57
58contains
59
67 function adjoint_output_init(precision, adjoint, adjoint_scalars, &
68 name, path) result(this)
69 integer, intent(inout) :: precision
70 class(adjoint_fluid_scheme_t), intent(in), target :: adjoint
71 class(adjoint_scalars_t), intent(in), optional, target :: adjoint_scalars
72 character(len=*), intent(in), optional :: name
73 character(len=*), intent(in), optional :: path
74 type(adjoint_output_t) :: this
75 character(len=1024) :: fname
76 integer :: i, n_scalars
77
78 if (present(name) .and. present(path)) then
79 fname = trim(path) // trim(name) // '.fld'
80 else if (present(name)) then
81 fname = trim(name) // '.fld'
82 else if (present(path)) then
83 fname = trim(path) // 'adjoint.fld'
84 else
85 fname = 'adjoint.fld'
86 end if
87
88 call this%init_base(fname, precision)
89
90 ! Calculate total number of fields
91 n_scalars = 0
92 if (present(adjoint_scalars)) then
93 n_scalars = size(adjoint_scalars%adjoint_scalar_fields)
94 end if
95
96 ! Initialize field list with appropriate size
97 call this%adjoint%init(4 + n_scalars)
98
99 call this%adjoint%assign(1, adjoint%p_adj)
100 call this%adjoint%assign(2, adjoint%u_adj)
101 call this%adjoint%assign(3, adjoint%v_adj)
102 call this%adjoint%assign(4, adjoint%w_adj)
103
104 ! Assign all scalar fields
105 if (present(adjoint_scalars)) then
106 do i = 1, n_scalars
107 call this%adjoint%assign(4 + i, &
108 adjoint_scalars%adjoint_scalar_fields(i)%s_adj)
109 end do
110 end if
111
112 end function adjoint_output_init
113
117 subroutine adjoint_output_sample(this, t)
118 class(adjoint_output_t), intent(inout) :: this
119 real(kind=rp), intent(in) :: t
120 integer :: i
121
122 if (neko_bcknd_device .eq. 1) then
123
124 associate(fields => this%adjoint%items)
125 do i = 1, size(fields)
126 call device_memcpy(fields(i)%ptr%x, fields(i)%ptr%x_d, &
127 fields(i)%ptr%dof%size(), device_to_host, &
128 sync = (i .eq. size(fields))) ! Sync on the last field
129 end do
130 end associate
131
132 end if
133
134 call this%file_%write(this%adjoint, t)
135
136 end subroutine adjoint_output_sample
137
138end 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.