Skip to content
GitLab
Projects
Groups
Snippets
/
Help
Help
Support
Community forum
Keyboard shortcuts
?
Submit feedback
Contribute to GitLab
Sign in
Toggle navigation
Menu
Open sidebar
SWIFT
SWIFTsim
Commits
d3340fdd
Commit
d3340fdd
authored
Apr 16, 2019
by
Matthieu Schaller
Browse files
Only do feedback for stars with an age > 0.
parent
fdb8eab7
Changes
7
Hide whitespace changes
Inline
Side-by-side
src/feedback/EAGLE/feedback.h
View file @
d3340fdd
...
...
@@ -56,9 +56,6 @@ __attribute__((always_inline)) INLINE static void feedback_init_spart(
__attribute__
((
always_inline
))
INLINE
static
void
feedback_first_init_spart
(
struct
spart
*
sp
,
const
struct
feedback_props
*
feedback_props
)
{
sp
->
birth_density
=
-
1
.
f
;
sp
->
birth_time
=
feedback_props
->
spart_first_init_birth_time
;
feedback_init_spart
(
sp
);
}
...
...
@@ -74,8 +71,31 @@ __attribute__((always_inline)) INLINE static void feedback_first_init_spart(
__attribute__
((
always_inline
))
INLINE
static
void
feedback_prepare_spart
(
struct
spart
*
sp
,
const
struct
feedback_props
*
feedback_props
)
{
/* Zero all the output */
bzero
(
&
sp
->
feedback_data
,
sizeof
(
struct
feedback_spart_data
));
/* Zero the amount of mass that is distributed */
sp
->
feedback_data
.
to_distribute
.
mass
=
0
.
f
;
/* Zero the number of SN */
sp
->
feedback_data
.
to_distribute
.
num_SNIa
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
num_SNII
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
num_SNe
=
0
.
f
;
/* Zero the enrichment quantities */
for
(
int
i
=
0
;
i
<
chemistry_element_count
;
i
++
)
sp
->
feedback_data
.
to_distribute
.
metal_mass
[
i
]
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
total_metal_mass
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
mass_from_AGB
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
metal_mass_from_AGB
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
mass_from_SNII
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
metal_mass_from_SNII
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
mass_from_SNIa
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
metal_mass_from_SNIa
=
0
.
f
;
sp
->
feedback_data
.
to_distribute
.
Fe_mass_from_SNIa
=
0
.
f
;
/* Zero the energy to inject */
sp
->
feedback_data
.
to_distribute
.
d_energy
=
0
.
f
;
/* Zero the feedback probability */
sp
->
feedback_data
.
to_distribute
.
heating_probability
=
0
.
f
;
}
/**
...
...
@@ -85,34 +105,22 @@ __attribute__((always_inline)) INLINE static void feedback_prepare_spart(
* this information to a different MPI rank.
*
* @param sp The particle to act upon
* @param feedback_propss The #feedback_props structure.
* @param cosmo The current cosmological model.
* @param stars_properties The #stars_props
* @param us The unit system.
* @param star_age_beg_step The age of the star at the star of the time-step in
* internal units.
* @param dt The time-step size of this star in internal units.
*/
__attribute__
((
always_inline
))
INLINE
static
void
feedback_evolve_spart
(
struct
spart
*
restrict
sp
,
const
struct
feedback_props
*
feedback_props
,
const
struct
cosmology
*
cosmo
,
const
struct
unit_system
*
us
,
double
star_age
,
double
dt
)
{
/* Zero the number of SN and amount of mass that is distributed */
sp
->
feedback_data
.
to_distribute
.
num_SNIa
=
0
;
sp
->
feedback_data
.
to_distribute
.
num_SNII
=
0
;
sp
->
feedback_data
.
to_distribute
.
mass
=
0
;
/* Zero the enrichment quantities */
for
(
int
i
=
0
;
i
<
chemistry_element_count
;
i
++
)
sp
->
feedback_data
.
to_distribute
.
metal_mass
[
i
]
=
0
;
sp
->
feedback_data
.
to_distribute
.
total_metal_mass
=
0
;
sp
->
feedback_data
.
to_distribute
.
mass_from_AGB
=
0
;
sp
->
feedback_data
.
to_distribute
.
metal_mass_from_AGB
=
0
;
sp
->
feedback_data
.
to_distribute
.
mass_from_SNII
=
0
;
sp
->
feedback_data
.
to_distribute
.
metal_mass_from_SNII
=
0
;
sp
->
feedback_data
.
to_distribute
.
mass_from_SNIa
=
0
;
sp
->
feedback_data
.
to_distribute
.
metal_mass_from_SNIa
=
0
;
sp
->
feedback_data
.
to_distribute
.
Fe_mass_from_SNIa
=
0
;
const
double
star_age_beg_step
,
const
double
dt
)
{
/* Compute amount of enrichment and feedback that needs to be done in this
* step */
compute_stellar_evolution
(
feedback_props
,
cosmo
,
sp
,
us
,
star_age
,
dt
);
compute_stellar_evolution
(
feedback_props
,
cosmo
,
sp
,
us
,
star_age_beg_step
,
dt
);
/* Decrease star mass by amount of mass distributed to gas neighbours */
sp
->
mass
-=
sp
->
feedback_data
.
to_distribute
.
mass
;
...
...
src/feedback/EAGLE/feedback_properties.h
View file @
d3340fdd
...
...
@@ -148,9 +148,6 @@ struct feedback_props {
/* wind delay time for SNII */
float
SNII_wind_delay
;
/* Value to set birth time of stars read from ICs */
float
spart_first_init_birth_time
;
};
void
feedback_props_init
(
struct
feedback_props
*
fp
,
...
...
src/runner.c
View file @
d3340fdd
...
...
@@ -205,19 +205,6 @@ void runner_do_stars_ghost(struct runner *r, struct cell *c, int timer) {
/* Get a direct pointer on the part. */
struct
spart
*
sp
=
&
sparts
[
sid
[
i
]];
const
integertime_t
ti_step
=
get_integer_timestep
(
sp
->
time_bin
);
const
integertime_t
ti_begin
=
get_integer_time_begin
(
e
->
ti_current
-
1
,
sp
->
time_bin
);
/* Get particle time-step */
double
dt
;
if
(
with_cosmology
)
{
dt
=
cosmology_get_delta_time
(
e
->
cosmology
,
ti_begin
,
ti_begin
+
ti_step
);
}
else
{
dt
=
get_timestep
(
sp
->
time_bin
,
e
->
time_base
);
}
#ifdef SWIFT_DEBUG_CHECKS
/* Is this part within the timestep? */
if
(
!
spart_is_active
(
sp
,
e
))
...
...
@@ -375,17 +362,42 @@ void runner_do_stars_ghost(struct runner *r, struct cell *c, int timer) {
/* Only do feedback if stars have a reasonable birth time */
if
(
sp
->
birth_time
!=
-
1
.)
{
/* Calculate age of the star */
double
star_age
;
const
integertime_t
ti_step
=
get_integer_timestep
(
sp
->
time_bin
);
const
integertime_t
ti_begin
=
get_integer_time_begin
(
e
->
ti_current
-
1
,
sp
->
time_bin
);
/* Get particle time-step */
double
dt
;
if
(
with_cosmology
)
{
star_age
=
cosmology_get_delta_time_from_scale_factors
(
dt
=
cosmology_get_delta_time
(
e
->
cosmology
,
ti_begin
,
ti_begin
+
ti_step
);
}
else
{
dt
=
get_timestep
(
sp
->
time_bin
,
e
->
time_base
);
}
/* Calculate age of the star at current time */
double
star_age_end_of_step
;
if
(
with_cosmology
)
{
star_age_end_of_step
=
cosmology_get_delta_time_from_scale_factors
(
cosmo
,
sp
->
birth_scale_factor
,
(
float
)
cosmo
->
a
);
}
else
{
star_age
=
e
->
time
-
sp
->
birth_time
;
star_age
_end_of_step
=
e
->
time
-
sp
->
birth_time
;
}
/* Compute the stellar evolution */
feedback_evolve_spart
(
sp
,
feedback_props
,
cosmo
,
us
,
star_age
,
dt
);
/* Reset the feedback fields of the star particle */
feedback_prepare_spart
(
sp
,
feedback_props
);
/* Has this star been around for a while ? */
if
(
star_age_end_of_step
>
0
.)
{
/* Age of the star at the start of the step */
const
double
star_age_beg_of_step
=
max
(
star_age_end_of_step
-
dt
,
0
.);
/* Compute the stellar evolution */
feedback_evolve_spart
(
sp
,
feedback_props
,
cosmo
,
us
,
star_age_beg_of_step
,
dt
);
}
}
}
...
...
src/stars/Default/stars_part.h
View file @
d3340fdd
...
...
@@ -87,7 +87,7 @@ struct spart {
/*! Feedback structure */
struct
feedback_spart_data
feedback_data
;
/*! Tracer structure */
struct
tracers_xpart_data
tracers_data
;
...
...
src/stars/EAGLE/stars.h
View file @
d3340fdd
...
...
@@ -66,6 +66,9 @@ __attribute__((always_inline)) INLINE static void stars_first_init_spart(
sp
->
time_bin
=
0
;
sp
->
birth_density
=
-
1
.
f
;
sp
->
birth_time
=
stars_properties
->
spart_first_init_birth_time
;
stars_init_spart
(
sp
);
}
...
...
src/stars/EAGLE/stars_io.h
View file @
d3340fdd
...
...
@@ -131,6 +131,11 @@ INLINE static void stars_props_init(struct stars_props *sp,
sp
->
log_max_h_change
=
p
->
log_max_h_change
;
else
sp
->
log_max_h_change
=
logf
(
powf
(
max_volume_change
,
hydro_dimension_inv
));
/* Read birth time to set all stars in ICs to (defaults to -1 to indicate star
* present in ICs) */
sp
->
spart_first_init_birth_time
=
parser_get_opt_param_float
(
params
,
"EAGLEFeedback:birth_time_override"
,
-
1
);
}
/**
...
...
src/stars/EAGLE/stars_part.h
View file @
d3340fdd
...
...
@@ -156,6 +156,9 @@ struct stars_props {
/*! Maximal change of h over one time-step */
float
log_max_h_change
;
/*! Value to set birth time of stars read from ICs */
float
spart_first_init_birth_time
;
};
#endif
/* SWIFT_EAGLE_STAR_PART_H */
Write
Preview
Supports
Markdown
0%
Try again
or
attach a new file
.
Cancel
You are about to add
0
people
to the discussion. Proceed with caution.
Finish editing this message first!
Cancel
Please
register
or
sign in
to comment