156 ,ig_hexa_gnode_xiloc&
157 ,ig_hexa_gnode_etloc&
158 ,ig_hexa_gnode_etloc&
159 ,ig_hexa_gnode_zeloc&
164 ,ig_hexa_gnode_glonum
170 integer,
intent( in) :: ihexa
171 real ,
intent( in) :: xi
172 real ,
intent( in) :: eta
173 real ,
intent( in) :: zeta
174 real ,
intent(out) :: x
175 real ,
intent(out) :: y
176 real ,
intent(out) :: z
191 do m = 1,ig_hexa_nnode
193 lag_xi =
lagrap_geom(ig_hexa_gnode_xiloc(m),xi ,ig_line_nnode)
194 lag_et =
lagrap_geom(ig_hexa_gnode_etloc(m),eta ,ig_line_nnode)
195 lag_ze =
lagrap_geom(ig_hexa_gnode_zeloc(m),zeta,ig_line_nnode)
197 xgnode = rg_gnode_x(ig_hexa_gnode_glonum(m,ihexa))
198 ygnode = rg_gnode_y(ig_hexa_gnode_glonum(m,ihexa))
199 zgnode = rg_gnode_z(ig_hexa_gnode_glonum(m,ihexa))
201 x = x + lag_xi*lag_et*lag_ze*xgnode
202 y = y + lag_xi*lag_et*lag_ze*ygnode
203 z = z + lag_xi*lag_et*lag_ze*zgnode
228 ,ig_quad_gnode_xiloc&
229 ,ig_quad_gnode_etloc&
234 ,ig_quadf_gnode_glonum
240 integer,
intent( in) :: iquad
241 real ,
intent( in) :: xi
242 real ,
intent( in) :: eta
243 real ,
intent(out) :: x
244 real ,
intent(out) :: y
245 real ,
intent(out) :: z
259 do m = 1,ig_quad_nnode
261 lag_xi =
lagrap_geom(ig_quad_gnode_xiloc(m),xi ,ig_line_nnode)
262 lag_et =
lagrap_geom(ig_quad_gnode_etloc(m),eta,ig_line_nnode)
264 xgnode = rg_gnode_x(ig_quadf_gnode_glonum(m,iquad))
265 ygnode = rg_gnode_y(ig_quadf_gnode_glonum(m,iquad))
266 zgnode = rg_gnode_z(ig_quadf_gnode_glonum(m,iquad))
268 x = x + lag_xi*lag_et*xgnode
269 y = y + lag_xi*lag_et*ygnode
270 z = z + lag_xi*lag_et*zgnode
293 ,ig_quadf_gnode_glonum&
295 ,ig_quad_gnode_xiloc&
302 integer,
intent(in) :: iquad
303 real ,
intent(in) :: xi
304 real ,
intent(in) :: eta
314 do m = 1,ig_quad_nnode
316 lag_xi =
lagrap_geom(ig_quad_gnode_xiloc(m),xi ,ig_line_nnode)
317 lag_et =
lagrap_geom(ig_quad_gnode_etloc(m),eta,ig_line_nnode)
319 zgnode = rg_gnode_z(ig_quadf_gnode_glonum(m,iquad))
321 z = z + lag_xi*lag_et*zgnode
This module defines all global variables of EFISPEC3D. Scalar variables are initialized directly in t...
real function, public lagrap_geom(i, x, n)
function to compute value of order Lagrange polynomial of geometric node i at abscissa : ...
This module contains subroutines to compute global -coordinates of a given point in the physical doma...
subroutine, public compute_hexa_point_coord(ihexa, xi, eta, zeta, x, y, z)
subroutine to compute -coordinates of a point knowing to which hexahedron element it belongs and its ...
real function, public compute_quad_point_coord_z(iquad, xi, eta)
function to compute -coordinates of a point knowing to which quadrangle element it belongs and its lo...
subroutine, public compute_quad_point_coord(iquad, xi, eta, x, y, z)
subroutine to compute -coordinates of a point knowing to which quadrangle element it belongs and its ...
This module contains functions to compute Lagrange polynomials and its derivatives; and to interpolat...