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
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
43 implicit none
44 private
45
47 type, public, extends(output_t) :: adjoint_output_t
48 type(field_list_t) :: adjoint
49 contains
50 procedure, pass(this) :: sample => adjoint_output_sample
51 end type adjoint_output_t
52
53 interface adjoint_output_t
54 module procedure adjoint_output_init
55 end interface adjoint_output_t
56
57contains
58
59 function adjoint_output_init(precision, adjoint, scalar, name, path) &
60 result(this)
61 integer, intent(inout) :: precision
62 class(adjoint_fluid_scheme_t), intent(in), target :: adjoint
63 class(adjoint_scalar_scheme_t), intent(in), optional, target :: scalar
64 character(len=*), intent(in), optional :: name
65 character(len=*), intent(in), optional :: path
66 type(adjoint_output_t) :: this
67 character(len=1024) :: fname
68
69 if (present(name) .and. present(path)) then
70 fname = trim(path) // trim(name) // '.fld'
71 else if (present(name)) then
72 fname = trim(name) // '.fld'
73 else if (present(path)) then
74 fname = trim(path) // 'adjoint.fld'
75 else
76 fname = 'adjoint.fld'
77 end if
78
79 call this%init_base(fname, precision)
80
81 if (present(scalar)) then
82 call this%adjoint%init(5)
83 else
84 call this%adjoint%init(4)
85 end if
86
87 call this%adjoint%assign(1, adjoint%p_adj)
88 call this%adjoint%assign(2, adjoint%u_adj)
89 call this%adjoint%assign(3, adjoint%v_adj)
90 call this%adjoint%assign(4, adjoint%w_adj)
91
92 if (present(scalar)) then
93 call this%adjoint%assign(5, scalar%s_adj)
94 end if
95
96 end function adjoint_output_init
97
99 subroutine adjoint_output_sample(this, t)
100 class(adjoint_output_t), intent(inout) :: this
101 real(kind=rp), intent(in) :: t
102 integer :: i
103
104 if (neko_bcknd_device .eq. 1) then
105
106 associate(fields => this%adjoint%items)
107 do i = 1, size(fields)
108 call device_memcpy(fields(i)%ptr%x, fields(i)%ptr%x_d, &
109 fields(i)%ptr%dof%size(), device_to_host, &
110 sync = (i .eq. size(fields))) ! Sync on the last field
111 end do
112 end associate
113
114 end if
115
116 call this%file_%write(this%adjoint, t)
117
118 end subroutine adjoint_output_sample
119
120end module adjoint_output
Defines an output for a adjoint.
Contains the adjoint_scalar_scheme_t type.
Base type for a scalar advection-diffusion solver.