diff --git a/generate_all_operators.py b/generate_all_operators.py index 2f81aeefce7ebb8c46278109dcb12d4b29bb5e3b..d1097c1008f7140c91007aef0d79e0f70db9fe94 100644 --- a/generate_all_operators.py +++ b/generate_all_operators.py @@ -43,6 +43,7 @@ from hog.forms import ( nonlinear_diffusion, nonlinear_diffusion_newton_galerkin, supg_diffusion, + grad_rho_by_rho_dot_u, ) from hog.forms_vectorial import curl_curl, curl_curl_plus_mass, mass_n1e1 from hog.function_space import ( @@ -615,6 +616,24 @@ def all_operators( ) ) + p2vec_grad_rho = partial( + grad_rho_by_rho_dot_u, + density_function_space=P2, + ) + ops.append( + OperatorInfo( + mapping=f"P2VectorToP1", + name=f"GradRhoByRhoDotU", + trial_space=P1, + test_space=P2Vector, + form=p2vec_grad_rho, + type_descriptor=type_descriptor, + geometries=list(geometries), + opts=opts, + blending=blending, + ) + ) + for c in [0, 1, 2]: # fmt: off if c == 2: diff --git a/hog/forms.py b/hog/forms.py index 8f316d6f6dfa1b085a763c5c7a805713775d8082..16b620ed41fd92f63a00c4ca33c9e9d8e0219d76 100644 --- a/hog/forms.py +++ b/hog/forms.py @@ -1818,6 +1818,85 @@ Weak formulation return Form(mat, tabulation, symmetric=False, docstring=docstring) +def grad_rho_by_rho_dot_u( + trial: FunctionSpace, + test: FunctionSpace, + geometry: ElementGeometry, + symbolizer: Symbolizer, + blending: GeometryMap = IdentityMap(), + density_function_space: Optional[FunctionSpace] = None, +) -> Form: + docstring = f""" +RHS operator for the frozen velocity approach. + +Geometry map: {blending} + +Weak formulation + + u: trial function (vectorial space: {trial}) + v: test function (space: {test}) + rho: coefficient (space: {density_function_space}) + + ∫ ((âˆ‡Ï / Ï) · u) v +""" + + with TimedLogger("assembling grad_rho_by_rho_dot_u-mass matrix", level=logging.DEBUG): + tabulation = Tabulation(symbolizer) + + jac_affine_det = symbolizer.abs_det_jac_ref_to_affine() + jac_affine_inv = symbolizer.jac_ref_to_affine_inv(geometry.dimensions) + + if isinstance(blending, ExternalMap): + jac_blending = jac_affine_to_physical(geometry, symbolizer) + else: + affine_coords = trafo_ref_to_affine(geometry, symbolizer) + # jac_blending = blending.jacobian(affine_coords) + jac_blending_inv = symbolizer.jac_affine_to_blending_inv(geometry.dimensions) + jac_blending_det = symbolizer.abs_det_jac_affine_to_blending() + + # jac_blending_det = abs(det(jac_blending)) + + mat = create_empty_element_matrix(trial, test, geometry) + + it = element_matrix_iterator(trial, test, geometry) + + if density_function_space: + rho, dof_symbols = fem_function_on_element( + density_function_space, + geometry, + symbolizer, + domain="reference", + function_id="rho", + ) + + grad_rho, _ = fem_function_gradient_on_element( + density_function_space, + geometry, + symbolizer, + domain="reference", + function_id="grad_rho", + dof_symbols=dof_symbols, + ) + + with TimedLogger( + f"integrating {mat.shape[0] * mat.shape[1]} expressions", + level=logging.DEBUG, + ): + for data in it: + phi = data.trial_shape + psi = data.test_shape + # print("phi = ", psi.shape) + # phi_psi_jac_affine_det = tabulation.register_factor( + # "phi_psi_jac_affine_det", + # sp.Matrix([phi * psi * jac_affine_det]), + # )[0] + if blending == IdentityMap(): + form = dot(((jac_affine_inv.T * grad_rho) / rho[0]), phi) * psi * jac_affine_det + else: + form = dot(((jac_blending_inv.T * jac_affine_inv.T * grad_rho) / rho[0]), phi) * psi * jac_affine_det * jac_blending_det + mat[data.row, data.col] = form + + return Form(mat, tabulation, symmetric=trial == test, docstring=docstring) def zero_form( trial: FunctionSpace, test: FunctionSpace, geometry: ElementGeometry