Neko-TOP
A portable framework for high-order spectral element flow toplogy optimization.
Loading...
Searching...
No Matches
heaviside_mapping.f90
Go to the documentation of this file.
1
34!
37 use num_types, only: rp
38 use mapping, only: mapping_t
39 use json_module, only: json_file
40 use field, only: field_t
41 use coefs, only: coef_t
42 use neko_config, only: neko_bcknd_device
43 use device_heaviside_mapping, only: &
44 device_heaviside_mapping_apply, &
45 device_heaviside_mapping_apply_backward
46 use heaviside_mapping_cpu, only: &
49 use json_utils, only: json_get_or_default
50 use utils, only: neko_error
51 use continuation_scheduler, only: nekotop_continuation
52
53 implicit none
54 private
55
65 type, public, extends(mapping_t) :: heaviside_mapping_t
67 real(kind=rp) :: beta
69 real(kind=rp) :: eta
70
71 contains
73 procedure, pass(this) :: init => heaviside_mapping_init_from_json
75 procedure, pass(this) :: init_from_attributes => &
76 heaviside_mapping_init_from_attributes
78 procedure, pass(this) :: free => heaviside_mapping_free
80 procedure, pass(this) :: forward_mapping => heaviside_mapping_forward
82 procedure, pass(this) :: backward_mapping => heaviside_mapping_backward
83 end type heaviside_mapping_t
84
85contains
86
88 subroutine heaviside_mapping_init_from_json(this, json, coef)
89 class(heaviside_mapping_t), intent(inout) :: this
90 type(json_file), intent(inout) :: json
91 type(coef_t), intent(inout) :: coef
92 real(kind=rp) :: eta, beta
93
94 call json_get_or_default(json, 'eta', eta, 0.5_rp)
95
96 call nekotop_continuation%json_get_or_register(json, 'beta', this%beta, &
97 beta)
98
99 ! Initialize base
100 call this%init_base(json, coef, "heaviside_mapping")
101 call this%init_from_attributes(coef, beta, eta)
103
105 subroutine heaviside_mapping_init_from_attributes(this, coef, beta, eta)
106 class(heaviside_mapping_t), intent(inout) :: this
107 type(coef_t), intent(inout) :: coef
108 real(kind=rp), intent(in) :: beta
109 real(kind=rp), intent(in) :: eta
110
111 if (beta .le. 0.0_rp) then
112 call neko_error('"beta" must be > 0 in heaviside mapping')
113 end if
114
115 if (eta .lt. 0.0_rp .or. eta .gt. 1.0_rp) then
116 call neko_error('"eta" must be in [0, 1] in heaviside mapping')
117 end if
118
119 this%beta = beta
120 this%eta = eta
121 end subroutine heaviside_mapping_init_from_attributes
122
124 subroutine heaviside_mapping_free(this)
125 class(heaviside_mapping_t), intent(inout) :: this
126
127 call this%free_base()
128 end subroutine heaviside_mapping_free
129
131 subroutine heaviside_mapping_forward(this, X_out, X_in)
132 class(heaviside_mapping_t), intent(inout) :: this
133 type(field_t), intent(in) :: X_in
134 type(field_t), intent(inout) :: X_out
135 call heaviside_mapping_apply(this%beta, this%eta, x_out, x_in)
136 end subroutine heaviside_mapping_forward
137
139 subroutine heaviside_mapping_backward(this, sens_out, sens_in, X_in)
140 class(heaviside_mapping_t), intent(inout) :: this
141 type(field_t), intent(in) :: X_in
142 type(field_t), intent(in) :: sens_in
143 type(field_t), intent(inout) :: sens_out
144 call heaviside_mapping_apply_backward(this%beta, this%eta, &
145 sens_out, sens_in, x_in)
146 end subroutine heaviside_mapping_backward
147
153 subroutine heaviside_mapping_apply(beta, eta, X_out, X_in)
154 real(kind=rp), intent(in) :: beta, eta
155 type(field_t), intent(in) :: x_in
156 type(field_t), intent(inout) :: X_out
157 integer :: n
158 real(kind=rp) :: den, beta_eta
159
160 beta_eta = beta * eta
161 den = tanh(beta_eta) + tanh(beta * (1.0_rp - eta))
162
163 if (abs(den) .le. tiny(den)) then
164 call neko_error('invalid denominator in heaviside mapping')
165 end if
166
167 n = x_in%dof%size()
168 if (neko_bcknd_device .eq. 1) then
169 call device_heaviside_mapping_apply(beta, eta, &
170 x_out%x_d, x_in%x_d, n)
171 else
172 call heaviside_mapping_apply_cpu(beta, eta, &
173 x_out%x, x_in%x, n)
174 end if
175 end subroutine heaviside_mapping_apply
176
183 subroutine heaviside_mapping_apply_backward(beta, eta, sens_out, sens_in, &
184 X_in)
185 real(kind=rp), intent(in) :: beta, eta
186 type(field_t), intent(in) :: x_in
187 type(field_t), intent(in) :: sens_in
188 type(field_t), intent(inout) :: sens_out
189 integer :: n
190 real(kind=rp) :: den, beta_eta
191
192 beta_eta = beta * eta
193 den = tanh(beta_eta) + tanh(beta * (1.0_rp - eta))
194
195 if (abs(den) .le. tiny(den)) then
196 call neko_error('invalid denominator in heaviside mapping')
197 end if
198
199 n = x_in%dof%size()
200 if (neko_bcknd_device .eq. 1) then
201 call device_heaviside_mapping_apply_backward(beta, eta, &
202 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
203 else
205 sens_out%x, sens_in%x, x_in%x, n)
206 end if
207 end subroutine heaviside_mapping_apply_backward
208
209end module heaviside_mapping
Continuation scheduler for the optimization loop.
CPU backend for smooth Heaviside mapping operations.
subroutine, public heaviside_mapping_apply_cpu(beta, eta, x_out, x_in, n)
Apply smooth Heaviside mapping on CPU.
subroutine, public heaviside_mapping_apply_backward_cpu(beta, eta, sens_out, sens_in, x_in, n)
Apply smooth Heaviside mapping chain rule on CPU.
Smooth Heaviside mapping.
subroutine heaviside_mapping_init_from_json(this, json, coef)
Constructor from json.
Mappings to be applied to a scalar field.
Definition mapping.f90:36
Base abstract class for mapping.
Definition mapping.f90:46