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, json_get_or_default
50 use utils, only: neko_error
51 implicit none
52 private
53
63 type, public, extends(mapping_t) :: heaviside_mapping_t
65 real(kind=rp) :: beta
67 real(kind=rp) :: eta
68
69 contains
71 procedure, pass(this) :: init => heaviside_mapping_init_from_json
73 procedure, pass(this) :: init_from_attributes => &
74 heaviside_mapping_init_from_attributes
76 procedure, pass(this) :: free => heaviside_mapping_free
78 procedure, pass(this) :: forward_mapping => heaviside_mapping_forward
80 procedure, pass(this) :: backward_mapping => heaviside_mapping_backward
81 end type heaviside_mapping_t
82
83contains
84
86 subroutine heaviside_mapping_init_from_json(this, json, coef)
87 class(heaviside_mapping_t), intent(inout) :: this
88 type(json_file), intent(inout) :: json
89 type(coef_t), intent(inout) :: coef
90 real(kind=rp) :: beta, eta
91
92 call json_get(json, 'beta', beta)
93 call json_get_or_default(json, 'eta', eta, 0.5_rp)
94
95 call this%init_base(json, coef)
96 call this%init_from_attributes(coef, beta, eta)
98
100 subroutine heaviside_mapping_init_from_attributes(this, coef, beta, eta)
101 class(heaviside_mapping_t), intent(inout) :: this
102 type(coef_t), intent(inout) :: coef
103 real(kind=rp), intent(in) :: beta
104 real(kind=rp), intent(in) :: eta
105
106 if (beta .le. 0.0_rp) then
107 call neko_error('"beta" must be > 0 in heaviside mapping')
108 end if
109
110 if (eta .lt. 0.0_rp .or. eta .gt. 1.0_rp) then
111 call neko_error('"eta" must be in [0, 1] in heaviside mapping')
112 end if
113
114 this%beta = beta
115 this%eta = eta
116 end subroutine heaviside_mapping_init_from_attributes
117
119 subroutine heaviside_mapping_free(this)
120 class(heaviside_mapping_t), intent(inout) :: this
121
122 call this%free_base()
123 end subroutine heaviside_mapping_free
124
126 subroutine heaviside_mapping_forward(this, X_out, X_in)
127 class(heaviside_mapping_t), intent(inout) :: this
128 type(field_t), intent(in) :: X_in
129 type(field_t), intent(inout) :: X_out
130 call heaviside_mapping_apply(this%beta, this%eta, x_out, x_in)
131 end subroutine heaviside_mapping_forward
132
134 subroutine heaviside_mapping_backward(this, sens_out, sens_in, X_in)
135 class(heaviside_mapping_t), intent(inout) :: this
136 type(field_t), intent(in) :: X_in
137 type(field_t), intent(in) :: sens_in
138 type(field_t), intent(inout) :: sens_out
139 call heaviside_mapping_apply_backward(this%beta, this%eta, &
140 sens_out, sens_in, x_in)
141 end subroutine heaviside_mapping_backward
142
148 subroutine heaviside_mapping_apply(beta, eta, X_out, X_in)
149 real(kind=rp), intent(in) :: beta, eta
150 type(field_t), intent(in) :: x_in
151 type(field_t), intent(inout) :: X_out
152 integer :: n
153 real(kind=rp) :: den, beta_eta
154
155 beta_eta = beta * eta
156 den = tanh(beta_eta) + tanh(beta * (1.0_rp - eta))
157
158 if (abs(den) .le. tiny(den)) then
159 call neko_error('invalid denominator in heaviside mapping')
160 end if
161
162 n = x_in%dof%size()
163 if (neko_bcknd_device .eq. 1) then
164 call device_heaviside_mapping_apply(beta, eta, &
165 x_out%x_d, x_in%x_d, n)
166 else
167 call heaviside_mapping_apply_cpu(beta, eta, &
168 x_out%x, x_in%x, n)
169 end if
170 end subroutine heaviside_mapping_apply
171
178 subroutine heaviside_mapping_apply_backward(beta, eta, sens_out, sens_in, &
179 X_in)
180 real(kind=rp), intent(in) :: beta, eta
181 type(field_t), intent(in) :: x_in
182 type(field_t), intent(in) :: sens_in
183 type(field_t), intent(inout) :: sens_out
184 integer :: n
185 real(kind=rp) :: den, beta_eta
186
187 beta_eta = beta * eta
188 den = tanh(beta_eta) + tanh(beta * (1.0_rp - eta))
189
190 if (abs(den) .le. tiny(den)) then
191 call neko_error('invalid denominator in heaviside mapping')
192 end if
193
194 n = x_in%dof%size()
195 if (neko_bcknd_device .eq. 1) then
196 call device_heaviside_mapping_apply_backward(beta, eta, &
197 sens_out%x_d, sens_in%x_d, x_in%x_d, n)
198 else
200 sens_out%x, sens_in%x, x_in%x, n)
201 end if
202 end subroutine heaviside_mapping_apply_backward
203
204end module heaviside_mapping
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