/*******************************************************************************
* This file is part of SWIFT.
* Copyright (c) 2012 Pedro Gonnet (pedro.gonnet@durham.ac.uk)
* Matthieu Schaller (schaller@strw.leidenuniv.nl)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Lesser General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with this program. If not, see .
*
******************************************************************************/
#ifndef SWIFT_GEAR_CHEMISTRY_IACT_H
#define SWIFT_GEAR_CHEMISTRY_IACT_H
/**
* @file GEAR/chemistry_iact.h
* @brief Smooth metal interaction functions following the GEAR version of
* smooth metalicity.
*
* The interactions computed here are the ones presented in Wiersma, Schaye et
* al. 2009
*/
/**
* @brief do chemistry computation after the runner_iact_density (symmetric
* version)
*
* @param r2 Comoving square distance between the two particles.
* @param dx Comoving vector separating both particles (pi - pj).
* @param hi Comoving smoothing-length of particle i.
* @param hj Comoving smoothing-length of particle j.
* @param pi First particle.
* @param pj Second particle.
* @param a Current scale factor.
* @param H Current Hubble parameter.
*/
__attribute__((always_inline)) INLINE static void runner_iact_chemistry(
const float r2, const float dx[3], const float hi, const float hj,
struct part *restrict pi, struct part *restrict pj, const float a,
const float H) {
struct chemistry_part_data *chi = &pi->chemistry_data;
struct chemistry_part_data *chj = &pj->chemistry_data;
float wi;
float wj;
/* Get r */
const float r = sqrtf(r2);
/* Compute the kernel function for pi */
const float ui = r / hi;
kernel_eval(ui, &wi);
/* Compute the kernel function for pj */
const float uj = r / hj;
kernel_eval(uj, &wj);
/* Compute contribution to the smooth metallicity */
for (int i = 0; i < GEAR_CHEMISTRY_ELEMENT_COUNT; i++) {
chi->smoothed_metal_mass_fraction[i] += chj->metal_mass[i] * wi;
chj->smoothed_metal_mass_fraction[i] += chi->metal_mass[i] * wj;
}
}
/**
* @brief do chemistry computation after the runner_iact_density (non symmetric
* version)
*
* @param r2 Comoving square distance between the two particles.
* @param dx Comoving vector separating both particles (pi - pj).
* @param hi Comoving smoothing-length of particle i.
* @param hj Comoving smoothing-length of particle j.
* @param pi First particle.
* @param pj Second particle (not updated).
* @param a Current scale factor.
* @param H Current Hubble parameter.
*/
__attribute__((always_inline)) INLINE static void runner_iact_nonsym_chemistry(
const float r2, const float dx[3], const float hi, const float hj,
struct part *restrict pi, const struct part *restrict pj, const float a,
const float H) {
struct chemistry_part_data *chi = &pi->chemistry_data;
const struct chemistry_part_data *chj = &pj->chemistry_data;
float wi;
/* Get r */
const float r = sqrtf(r2);
/* Compute the kernel function for pi */
const float ui = r / hi;
kernel_eval(ui, &wi);
/* Compute contribution to the smooth metallicity */
for (int i = 0; i < GEAR_CHEMISTRY_ELEMENT_COUNT; i++) {
chi->smoothed_metal_mass_fraction[i] += chj->metal_mass[i] * wi;
}
}
/**
* @brief do metal diffusion computation in the force loop
* (symmetric version)
*
* @param r2 Comoving square distance between the two particles.
* @param dx Comoving vector separating both particles (pi - pj).
* @param hi Comoving smoothing-length of particle i.
* @param hj Comoving smoothing-length of particle j.
* @param pi First particle.
* @param pj Second particle.
* @param a Current scale factor.
* @param H Current Hubble parameter.
*/
__attribute__((always_inline)) INLINE static void
runner_iact_gradient_diffusion(const float r2, const float dx[3],
const float hi, const float hj,
struct part *restrict pi,
struct part *restrict pj, const float a,
const float H) {}
/**
* @brief do metal diffusion computation in the
* (nonsymmetric version)
*
* @param r2 Comoving square distance between the two particles.
* @param dx Comoving vector separating both particles (pi - pj).
* @param hi Comoving smoothing-length of particle i.
* @param hj Comoving smoothing-length of particle j.
* @param pi First particle.
* @param pj Second particle.
* @param a Current scale factor.
* @param H Current Hubble parameter.
*/
__attribute__((always_inline)) INLINE static void
runner_iact_nonsym_gradient_diffusion(const float r2, const float dx[3],
const float hi, const float hj,
struct part *restrict pi,
struct part *restrict pj, const float a,
const float H) {}
/**
* @brief Do metal diffusion computation in the (symmetric version)
*
* @param r2 Comoving square distance between the two particles.
* @param dx Comoving vector separating both particles (pi - pj).
* @param hi Comoving smoothing-length of particle i.
* @param hj Comoving smoothing-length of particle j.
* @param pi First particle.
* @param pj Second particle.
* @param a Current scale factor.
* @param H Current Hubble parameter.
* @param time_base The time base used in order to convert integer to float
* time.
* @param ti_current The current time (in integer)
* @param cosmo The #cosmology.
* @param with_cosmology Are we running with cosmology?
* @param chem_data The global properties of the chemistry scheme.
*
*/
__attribute__((always_inline)) INLINE static void runner_iact_diffusion(
const float r2, const float dx[3], const float hi, const float hj,
struct part *restrict pi, struct part *restrict pj, const float a,
const float H, const float time_base, const integertime_t t_current,
const struct cosmology *cosmo, const int with_cosmology,
const struct chemistry_global_data *chem_data) {}
/**
* @brief do metal diffusion computation in the force loop
* (nonsymmetric version)
*
* @param r2 Comoving square distance between the two particles.
* @param dx Comoving vector separating both particles (pi - pj).
* @param hi Comoving smoothing-length of particle i.
* @param hj Comoving smoothing-length of particle j.
* @param pi First particle.
* @param pj Second particle.
* @param a Current scale factor.
* @param H Current Hubble parameter.
* @param time_base The time base used in order to convert integer to float
* time.
* @param ti_current The current time (in integer)
* @param cosmo The #cosmology.
* @param with_cosmology Are we running with cosmology?
* @param chem_data The global properties of the chemistry scheme.
*
*/
__attribute__((always_inline)) INLINE static void runner_iact_nonsym_diffusion(
const float r2, const float dx[3], const float hi, const float hj,
struct part *restrict pi, struct part *restrict pj, const float a,
const float H, const float time_base, const integertime_t t_current,
const struct cosmology *cosmo, const int with_cosmology,
const struct chemistry_global_data *chem_data) {}
#endif /* SWIFT_GEAR_CHEMISTRY_IACT_H */