Page Menu
Home
c4science
Search
Configure Global Search
Log In
Files
F90824354
chimie.c
No One
Temporary
Actions
Download File
Edit File
Delete File
View Transforms
Subscribe
Mute Notifications
Award Token
Subscribers
None
File Metadata
Details
File Info
Storage
Attached
Created
Tue, Nov 5, 02:50
Size
62 KB
Mime Type
text/x-c
Expires
Thu, Nov 7, 02:50 (2 d)
Engine
blob
Format
Raw Data
Handle
22141080
Attached To
rPNBODY pNbody
chimie.c
View Options
This document is not UTF8. It was detected as ISO-8859-1 (Latin 1) and converted to UTF8 for display.
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <math.h>
#include <mpi.h>
#include <gsl/gsl_math.h>
#include "allvars.h"
#include "proto.h"
#ifdef CHIMIE
/*! \file hydra.c
* \brief Computation of SPH forces and rate of entropy generation
*
* This file contains the "second SPH loop", where the SPH forces are
* computed, and where the rate of change of entropy due to the shock heating
* (via artificial viscosity) is computed.
*/
static
double
hubble_a
,
atime
,
hubble_a2
,
fac_mu
,
fac_vsic_fix
,
a3inv
,
fac_egy
;
#ifdef FEEDBACK
static
double
fac_pow
;
#endif
#ifdef PERIODIC
static
double
boxSize
,
boxHalf
;
#ifdef LONG_X
static
double
boxSize_X
,
boxHalf_X
;
#else
#define boxSize_X boxSize
#define boxHalf_X boxHalf
#endif
#ifdef LONG_Y
static
double
boxSize_Y
,
boxHalf_Y
;
#else
#define boxSize_Y boxSize
#define boxHalf_Y boxHalf
#endif
#ifdef LONG_Z
static
double
boxSize_Z
,
boxHalf_Z
;
#else
#define boxSize_Z boxSize
#define boxHalf_Z boxHalf
#endif
#endif
/****************************************************************************************/
/*
/*
/*
/* GADGET CHIMIE PART
/*
/*
/*
/****************************************************************************************/
#define MAXPTS 10
#define MAXDATASIZE 200
#define KPC_IN_CM 3.085e+21
static
int
verbose
=
0
;
static
double
*
MassFracSNII
;
static
double
*
SingleMassFracSNII
;
static
double
*
EjectedMass
;
static
double
*
SingleEjectedMass
;
static
double
**
MassFracSNIIs
;
static
double
**
SingleMassFracSNIIs
;
static
double
**
EjectedMasss
;
static
double
**
SingleEjectedMasss
;
/* intern global variables */
static
struct
local_params_chimie
{
float
coeff_z
[
3
][
3
];
float
Mmin
,
Mmax
;
int
n
;
float
ms
[
MAXPTS
];
float
as
[
MAXPTS
+
1
];
float
bs
[
MAXPTS
+
1
];
float
fs
[
MAXPTS
];
double
imf_Ntot
;
float
SNII_Mmin
;
float
SNII_Mmax
;
float
SNII_cte
;
float
SNII_a
;
float
SNIa_Mpl
;
float
SNIa_Mpu
;
float
SNIa_a
;
float
SNIa_cte
;
float
SNIa_Mdl1
;
float
SNIa_Mdu1
;
float
SNIa_a1
;
float
SNIa_b1
;
float
SNIa_cte1
;
float
SNIa_bb1
;
float
SNIa_Mdl2
;
float
SNIa_Mdu2
;
float
SNIa_a2
;
float
SNIa_b2
;
float
SNIa_cte2
;
float
SNIa_bb2
;
float
Mco
;
int
npts
;
int
nelts
;
}
*
Cps
,
*
Cp
;
static
struct
local_elts_chimie
{
float
Mmin
;
/* minimal mass */
float
Step
;
/* log of mass step */
float
Array
[
MAXDATASIZE
];
/* data */
float
Metal
[
MAXDATASIZE
];
/* data */
float
MSNIa
;
float
SolarAbundance
;
char
label
[
72
];
}
**
Elts
,
*
Elt
;
void
allocate_chimie
()
{
int
j
;
/* allocate Cp */
Cps
=
malloc
((
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
struct
local_params_chimie
));
/* allocate elts */
Elts
=
malloc
((
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
struct
local_elts_chimie
));
//for (j=0;j<All.ChimieNumberOfParameterFiles;j++)
// Elt[j] = malloc((nelts) * sizeof(struct local_elts_chimie));
MassFracSNIIs
=
malloc
((
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
double
));
EjectedMasss
=
malloc
((
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
double
));
SingleMassFracSNIIs
=
malloc
((
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
double
));
SingleEjectedMasss
=
malloc
((
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
double
));
}
void
allocate_Elts
(
int
it
)
{
/* allocate memory for elts */
if
(
Cps
[
it
].
npts
<=
MAXDATASIZE
)
{
Elts
[
it
]
=
malloc
((
Cps
[
it
].
nelts
+
2
)
*
sizeof
(
struct
local_elts_chimie
));
}
else
{
printf
(
"
\n
Cps[it].npts = %d > MAXDATASIZE = %d !!!
\n\n
"
,
Cps
[
it
].
npts
,
MAXDATASIZE
);
endrun
(
88800
);
}
}
void
set_table
(
int
i
)
{
if
(
i
>=
All
.
ChimieNumberOfParameterFiles
)
{
printf
(
"
\n
set_table : i>= %d !!!
\n\n
"
,
All
.
ChimieNumberOfParameterFiles
);
endrun
(
88809
);
}
else
{
Cp
=
&
Cps
[
i
];
Elt
=
Elts
[
i
];
MassFracSNII
=
MassFracSNIIs
[
i
];
SingleMassFracSNII
=
SingleMassFracSNIIs
[
i
];
EjectedMass
=
EjectedMasss
[
i
];
SingleEjectedMass
=
SingleEjectedMasss
[
i
];
}
}
void
read_chimie
(
char
*
filename
,
int
it
)
{
char
line
[
72
];
FILE
*
fd
;
int
i
,
j
;
if
(
ThisTask
==
0
)
{
printf
(
"reading %s ...
\n
"
,
filename
);
fd
=
fopen
(
filename
,
"r"
);
/* read Lifetime */
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g %g %g
\n
"
,
&
Cps
[
it
].
coeff_z
[
0
][
0
],
&
Cps
[
it
].
coeff_z
[
0
][
1
],
&
Cps
[
it
].
coeff_z
[
0
][
2
]);
fscanf
(
fd
,
"%g %g %g
\n
"
,
&
Cps
[
it
].
coeff_z
[
1
][
0
],
&
Cps
[
it
].
coeff_z
[
1
][
1
],
&
Cps
[
it
].
coeff_z
[
1
][
2
]);
fscanf
(
fd
,
"%g %g %g
\n
"
,
&
Cps
[
it
].
coeff_z
[
2
][
0
],
&
Cps
[
it
].
coeff_z
[
2
][
1
],
&
Cps
[
it
].
coeff_z
[
2
][
2
]);
fgets
(
line
,
sizeof
(
line
),
fd
);
/* IMF Parameters */
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g %g
\n
"
,
&
Cps
[
it
].
Mmin
,
&
Cps
[
it
].
Mmax
);
fscanf
(
fd
,
"%d
\n
"
,
&
Cps
[
it
].
n
);
if
(
Cps
[
it
].
n
>
0
)
for
(
i
=
0
;
i
<
Cps
[
it
].
n
;
i
++
)
fscanf
(
fd
,
"%g"
,
&
Cps
[
it
].
ms
[
i
]);
else
fgets
(
line
,
sizeof
(
line
),
fd
);
for
(
i
=
0
;
i
<
Cps
[
it
].
n
+
1
;
i
++
)
fscanf
(
fd
,
"%g"
,
&
Cps
[
it
].
as
[
i
]);
fgets
(
line
,
sizeof
(
line
),
fd
);
/* Parameters for SNII Rates */
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g
\n
"
,
&
Cps
[
it
].
SNII_Mmin
);
fgets
(
line
,
sizeof
(
line
),
fd
);
/* Parameters for SNIa Rates */
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g %g
\n
"
,
&
Cps
[
it
].
SNIa_Mpl
,
&
Cps
[
it
].
SNIa_Mpu
);
fscanf
(
fd
,
"%g
\n
"
,
&
Cps
[
it
].
SNIa_a
);
fscanf
(
fd
,
"%g %g %g
\n
"
,
&
Cps
[
it
].
SNIa_Mdl1
,
&
Cps
[
it
].
SNIa_Mdu1
,
&
Cps
[
it
].
SNIa_bb1
);
fscanf
(
fd
,
"%g %g %g
\n
"
,
&
Cps
[
it
].
SNIa_Mdl2
,
&
Cps
[
it
].
SNIa_Mdu2
,
&
Cps
[
it
].
SNIa_bb2
);
fgets
(
line
,
sizeof
(
line
),
fd
);
/* Metal injection SNII */
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%d %d
\n
"
,
&
Cps
[
it
].
npts
,
&
Cps
[
it
].
nelts
);
/* allocate memory for elts */
allocate_Elts
(
it
);
/* injected metals */
for
(
i
=
0
;
i
<
Cps
[
it
].
nelts
+
2
;
i
++
)
{
fgets
(
line
,
sizeof
(
line
),
fd
);
/* strip trailing line */
for
(
j
=
0
;
j
<
strlen
(
line
);
j
++
)
if
(
line
[
j
]
==
'\n'
||
line
[
j
]
==
'\r'
)
line
[
j
]
=
'\0'
;
/* copy labels */
strcpy
(
Elts
[
it
][
i
].
label
,
line
);
strcpy
(
Elts
[
it
][
i
].
label
,
&
Elts
[
it
][
i
].
label
[
2
]);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g %g
\n
"
,
&
Elts
[
it
][
i
].
Mmin
,
&
Elts
[
it
][
i
].
Step
);
for
(
j
=
0
;
j
<
Cps
[
it
].
npts
;
j
++
)
{
fscanf
(
fd
,
"%g
\n
"
,
&
Elts
[
it
][
i
].
Metal
[
j
]);
}
}
/* integrals of injected metals */
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%d %d
\n
"
,
&
Cps
[
it
].
npts
,
&
Cps
[
it
].
nelts
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
/* integrals of injected metals */
for
(
i
=
0
;
i
<
Cps
[
it
].
nelts
+
2
;
i
++
)
{
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g %g
\n
"
,
&
Elts
[
it
][
i
].
Mmin
,
&
Elts
[
it
][
i
].
Step
);
for
(
j
=
0
;
j
<
Cps
[
it
].
npts
;
j
++
)
{
fscanf
(
fd
,
"%g
\n
"
,
&
Elts
[
it
][
i
].
Array
[
j
]);
}
}
/* Metal injection SNIa */
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%g
\n
"
,
&
Cps
[
it
].
Mco
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
int
nelts
;
char
label
[
72
];
fscanf
(
fd
,
"%d
\n
"
,
&
nelts
);
/* check */
if
(
nelts
!=
Cps
[
it
].
nelts
)
{
printf
(
"
\n
The number of elements in SNII (=%d) is not identical to the on of SNIa (=%d) !!!
\n\n
"
,
Cps
[
it
].
nelts
,
nelts
);
printf
(
"This is not supported by the current implementation !!!
\n
"
);
endrun
(
88805
);
}
for
(
i
=
0
;
i
<
Cps
[
it
].
nelts
+
2
;
i
++
)
{
fgets
(
line
,
sizeof
(
line
),
fd
);
/* label */
/* check label */
/* strip trailing line */
for
(
j
=
0
;
j
<
strlen
(
line
);
j
++
)
if
(
line
[
j
]
==
'\n'
||
line
[
j
]
==
'\r'
)
line
[
j
]
=
'\0'
;
strcpy
(
label
,
line
);
strcpy
(
label
,
&
label
[
2
]);
if
(
strcmp
(
label
,
Elts
[
it
][
i
].
label
)
!=
0
)
{
printf
(
"
\n
Label of SNII element %d (=%s) is different from the SNIa one (=%s) !!!
\n\n
"
,
i
,
Elts
[
it
][
i
].
label
,
label
);
endrun
(
88806
);
}
//fgets(line, sizeof(line), fd);
fscanf
(
fd
,
"%g
\n
"
,
&
Elts
[
it
][
i
].
MSNIa
);
}
/* Solar Abundances */
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fgets
(
line
,
sizeof
(
line
),
fd
);
fscanf
(
fd
,
"%d
\n
"
,
&
nelts
);
/* check */
if
(
nelts
!=
Cps
[
it
].
nelts
)
{
printf
(
"
\n
The number of elements in SolarAbundances (=%d) is not identical to the on of SNIa (=%d) !!!
\n\n
"
,
Cps
[
it
].
nelts
,
nelts
);
printf
(
"This is not supported by the current implementation !!!
\n
"
);
endrun
(
88805
);
}
for
(
i
=
0
;
i
<
Cps
[
it
].
nelts
;
i
++
)
{
fgets
(
line
,
sizeof
(
line
),
fd
);
/* label */
/* check label */
/* strip trailing line */
for
(
j
=
0
;
j
<
strlen
(
line
);
j
++
)
if
(
line
[
j
]
==
'\n'
||
line
[
j
]
==
'\r'
)
line
[
j
]
=
'\0'
;
strcpy
(
label
,
line
);
strcpy
(
label
,
&
label
[
2
]);
if
(
strcmp
(
label
,
Elts
[
it
][
i
+
2
].
label
)
!=
0
)
{
printf
(
"
\n
Label of SNII element %d (=%s) is different from the SNIa one (=%s) !!!
\n\n
"
,
i
,
Elts
[
it
][
i
+
2
].
label
,
label
);
endrun
(
88806
);
}
//fgets(line, sizeof(line), fd);
fscanf
(
fd
,
"%g
\n
"
,
&
Elts
[
it
][
i
+
2
].
SolarAbundance
);
}
fclose
(
fd
);
}
/* send Cps */
MPI_Bcast
(
Cps
,
(
All
.
ChimieNumberOfParameterFiles
)
*
sizeof
(
struct
local_params_chimie
),
MPI_BYTE
,
0
,
MPI_COMM_WORLD
);
/* slaves allocate Elts */
if
(
ThisTask
!=
0
)
allocate_Elts
(
it
);
MPI_Bcast
(
Elts
[
it
],
(
Cps
[
it
].
nelts
+
2
)
*
sizeof
(
struct
local_elts_chimie
),
MPI_BYTE
,
0
,
MPI_COMM_WORLD
);
/* allocate memory */
MassFracSNIIs
[
it
]
=
malloc
((
Cps
[
it
].
nelts
+
2
)
*
sizeof
(
double
));
EjectedMasss
[
it
]
=
malloc
((
Cps
[
it
].
nelts
+
2
)
*
sizeof
(
double
));
SingleMassFracSNIIs
[
it
]
=
malloc
((
Cps
[
it
].
nelts
+
2
)
*
sizeof
(
double
));
SingleEjectedMasss
[
it
]
=
malloc
((
Cps
[
it
].
nelts
+
2
)
*
sizeof
(
double
));
if
(
verbose
&&
ThisTask
==
0
)
{
printf
(
"%g %g %g
\n
"
,
Cps
[
it
].
coeff_z
[
0
][
0
],
Cps
[
it
].
coeff_z
[
0
][
1
],
Cps
[
it
].
coeff_z
[
0
][
2
]);
printf
(
"%g %g %g
\n
"
,
Cps
[
it
].
coeff_z
[
1
][
0
],
Cps
[
it
].
coeff_z
[
1
][
1
],
Cps
[
it
].
coeff_z
[
1
][
2
]);
printf
(
"%g %g %g
\n
"
,
Cps
[
it
].
coeff_z
[
2
][
0
],
Cps
[
it
].
coeff_z
[
2
][
1
],
Cps
[
it
].
coeff_z
[
2
][
2
]);
printf
(
"
\n
"
);
printf
(
"
\n
IMF
\n
"
);
printf
(
"%g %g
\n
"
,
Cps
[
it
].
Mmin
,
Cps
[
it
].
Mmax
);
printf
(
"%d
\n
"
,
Cps
[
it
].
n
);
for
(
i
=
0
;
i
<
Cps
[
it
].
n
;
i
++
)
printf
(
"ms : %g "
,
Cps
[
it
].
ms
[
i
]);
printf
(
"
\n
"
);
for
(
i
=
0
;
i
<
Cps
[
it
].
n
+
1
;
i
++
)
printf
(
"as : %g "
,
Cps
[
it
].
as
[
i
]);
printf
(
"
\n
"
);
printf
(
"
\n
Rate SNII
\n
"
);
printf
(
"%g "
,
Cps
[
it
].
SNII_Mmin
);
printf
(
"
\n
"
);
printf
(
"
\n
Rate SNIa
\n
"
);
printf
(
"%g %g
\n
"
,
Cps
[
it
].
SNIa_Mpl
,
Cps
[
it
].
SNIa_Mpu
);
printf
(
"%g
\n
"
,
Cps
[
it
].
SNIa_a
);
printf
(
"%g %g %g
\n
"
,
Cps
[
it
].
SNIa_Mdl1
,
Cps
[
it
].
SNIa_Mdu1
,
Cps
[
it
].
SNIa_b1
);
printf
(
"%g %g %g
\n
"
,
Cps
[
it
].
SNIa_Mdl2
,
Cps
[
it
].
SNIa_Mdu2
,
Cps
[
it
].
SNIa_b2
);
printf
(
"
\n
"
);
for
(
i
=
0
;
i
<
Cps
[
it
].
nelts
+
2
;
i
++
)
{
printf
(
"> %g %g
\n
"
,
Elts
[
it
][
i
].
Mmin
,
Elts
[
it
][
i
].
Step
);
for
(
j
=
0
;
j
<
Cps
[
it
].
npts
;
j
++
)
{
printf
(
" %g
\n
"
,
Elts
[
it
][
i
].
Array
[
j
]);
}
}
printf
(
"
\n
"
);
printf
(
"%g
\n
"
,
Cps
[
it
].
Mco
);
for
(
i
=
0
;
i
<
Cps
[
it
].
nelts
+
2
;
i
++
)
printf
(
"%g
\n
"
,
Elts
[
it
][
i
].
MSNIa
);
printf
(
"
\n
"
);
}
}
/*
This function returns the mass fraction of a star of mass m
using the current IMF
*/
static
double
get_imf
(
double
m
)
{
int
i
;
int
n
;
n
=
Cp
->
n
;
/* convert m in msol */
m
=
m
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
if
(
n
==
0
)
return
Cp
->
bs
[
0
]
*
pow
(
m
,
Cp
->
as
[
0
]);
else
{
for
(
i
=
0
;
i
<
n
;
i
++
)
if
(
m
<
Cp
->
ms
[
i
])
return
Cp
->
bs
[
i
]
*
pow
(
m
,
Cp
->
as
[
i
]);
return
Cp
->
bs
[
n
]
*
pow
(
m
,
Cp
->
as
[
n
]);
}
}
/*
This function returns the mass fraction between m1 and m2
per mass unit, using the current IMF
*/
static
double
get_imf_M
(
double
m1
,
double
m2
)
{
int
i
;
int
n
;
double
p
;
double
integral
=
0
;
double
mmin
,
mmax
;
n
=
Cp
->
n
;
/* convert m in msol */
m1
=
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
m2
=
m2
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
if
(
n
==
0
)
{
p
=
Cp
->
as
[
0
]
+
1
;
integral
=
(
Cp
->
bs
[
0
]
/
p
)
*
(
pow
(
m2
,
p
)
-
pow
(
m1
,
p
)
);
//printf("--> %g %g %g %g int=%g\n",m1,m2,pow(m2,p), pow(m1,p),integral);
}
else
{
integral
=
0
;
/* first */
if
(
m1
<
Cp
->
ms
[
0
])
{
mmin
=
m1
;
mmax
=
dmin
(
Cp
->
ms
[
0
],
m2
);
p
=
Cp
->
as
[
0
]
+
1
;
integral
+=
(
Cp
->
bs
[
0
]
/
p
)
*
(
pow
(
mmax
,
p
)
-
pow
(
mmin
,
p
)
);
}
/* last */
if
(
m2
>
Cp
->
ms
[
n
-
1
])
{
mmin
=
dmax
(
Cp
->
ms
[
n
-
1
],
m1
);
mmax
=
m2
;
p
=
Cp
->
as
[
n
]
+
1
;
integral
+=
(
Cp
->
bs
[
n
]
/
p
)
*
(
pow
(
mmax
,
p
)
-
pow
(
mmin
,
p
)
);
}
/* loop over other segments */
for
(
i
=
0
;
i
<
n
-
1
;
i
++
)
{
mmin
=
dmax
(
Cp
->
ms
[
i
],
m1
);
mmax
=
dmin
(
Cp
->
ms
[
i
+
1
],
m2
);
if
(
mmin
<
mmax
)
{
p
=
Cp
->
as
[
i
+
1
]
+
1
;
integral
+=
(
Cp
->
bs
[
i
+
1
]
/
p
)
*
(
pow
(
mmax
,
p
)
-
pow
(
mmin
,
p
)
);
}
}
}
/* convert into mass unit mass unit */
/* integral = integral * SOLAR_MASS/All.UnitMass_in_g;*/
return
integral
;
}
/*
This function returns the number fraction between m1 and m2
per mass unit, using the current IMF
*/
static
double
get_imf_N
(
double
m1
,
double
m2
)
{
int
i
;
int
n
;
double
p
;
double
integral
=
0
;
double
mmin
,
mmax
;
n
=
Cp
->
n
;
/* convert m in msol */
m1
=
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
m2
=
m2
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
if
(
n
==
0
)
{
p
=
Cp
->
as
[
0
];
integral
=
(
Cp
->
bs
[
0
]
/
p
)
*
(
pow
(
m2
,
p
)
-
pow
(
m1
,
p
)
);
}
else
{
integral
=
0
;
/* first */
if
(
m1
<
Cp
->
ms
[
0
])
{
mmin
=
m1
;
mmax
=
dmin
(
Cp
->
ms
[
0
],
m2
);
p
=
Cp
->
as
[
0
];
integral
+=
(
Cp
->
bs
[
0
]
/
p
)
*
(
pow
(
mmax
,
p
)
-
pow
(
mmin
,
p
)
);
}
/* last */
if
(
m2
>
Cp
->
ms
[
n
-
1
])
{
mmin
=
dmax
(
Cp
->
ms
[
n
-
1
],
m1
);
mmax
=
m2
;
p
=
Cp
->
as
[
n
];
integral
+=
(
Cp
->
bs
[
n
]
/
p
)
*
(
pow
(
mmax
,
p
)
-
pow
(
mmin
,
p
)
);
}
/* loop over other segments */
for
(
i
=
0
;
i
<
n
-
1
;
i
++
)
{
mmin
=
dmax
(
Cp
->
ms
[
i
],
m1
);
mmax
=
dmin
(
Cp
->
ms
[
i
+
1
],
m2
);
if
(
mmin
<
mmax
)
{
p
=
Cp
->
as
[
i
+
1
];
integral
+=
(
Cp
->
bs
[
i
+
1
]
/
p
)
*
(
pow
(
mmax
,
p
)
-
pow
(
mmin
,
p
)
);
}
}
}
/* convert into mass unit mass unit */
integral
=
integral
/
SOLAR_MASS
*
All
.
UnitMass_in_g
;
return
integral
;
}
/*
This function returns the number fraction between m1 and m2
per mass unit, using the current IMF
*/
static
double
imf_sampling
()
{
int
i
;
int
n
;
double
m
;
double
f
;
double
pmin
,
pmax
;
n
=
Cp
->
n
;
/* init random */
//srandom(irand);
f
=
(
double
)
random
()
/
(
double
)
RAND_MAX
;
if
(
n
==
0
)
{
pmin
=
pow
(
Cp
->
Mmin
,
Cp
->
as
[
0
]);
pmax
=
pow
(
Cp
->
Mmax
,
Cp
->
as
[
0
]);
m
=
pow
(
f
*
(
pmax
-
pmin
)
+
pmin
,
1.
/
Cp
->
as
[
0
]);
return
m
*
SOLAR_MASS
/
All
.
UnitMass_in_g
;
}
else
{
if
(
f
<
Cp
->
fs
[
0
])
{
pmin
=
pow
(
Cp
->
Mmin
,
Cp
->
as
[
0
]);
m
=
pow
(
Cp
->
imf_Ntot
*
Cp
->
as
[
0
]
/
Cp
->
bs
[
0
]
*
(
f
-
0
)
+
pmin
,
1.
/
Cp
->
as
[
0
]);
return
m
*
SOLAR_MASS
/
All
.
UnitMass_in_g
;
}
for
(
i
=
0
;
i
<
n
-
1
;
i
++
)
{
if
(
f
<
Cp
->
fs
[
i
+
1
])
{
pmin
=
pow
(
Cp
->
ms
[
i
]
,
Cp
->
as
[
i
+
1
]);
m
=
pow
(
Cp
->
imf_Ntot
*
Cp
->
as
[
i
+
1
]
/
Cp
->
bs
[
i
+
1
]
*
(
f
-
Cp
->
fs
[
i
])
+
pmin
,
1.
/
Cp
->
as
[
i
+
1
]);
return
m
*
SOLAR_MASS
/
All
.
UnitMass_in_g
;
}
}
/* last portion */
pmin
=
pow
(
Cp
->
ms
[
n
-
1
]
,
Cp
->
as
[
n
]);
m
=
pow
(
Cp
->
imf_Ntot
*
Cp
->
as
[
n
]
/
Cp
->
bs
[
n
]
*
(
f
-
Cp
->
fs
[
n
-
1
])
+
pmin
,
1.
/
Cp
->
as
[
n
]);
return
m
*
SOLAR_MASS
/
All
.
UnitMass_in_g
;
}
}
/*
This function initialized the imf parameters
defined in the chimie file
*/
void
init_imf
(
void
)
{
float
integral
=
0
;
float
p
;
float
cte
;
int
i
,
n
;
double
mmin
,
mmax
;
n
=
Cp
->
n
;
if
(
n
==
0
)
{
p
=
Cp
->
as
[
0
]
+
1
;
integral
=
integral
+
(
pow
(
Cp
->
Mmax
,
p
)
-
pow
(
Cp
->
Mmin
,
p
))
/
(
p
)
;
Cp
->
bs
[
0
]
=
1.
/
integral
;
}
else
{
cte
=
1.0
;
if
(
Cp
->
Mmin
<
Cp
->
ms
[
0
])
{
p
=
Cp
->
as
[
0
]
+
1
;
integral
=
integral
+
(
pow
(
Cp
->
ms
[
0
],
p
)
-
pow
(
Cp
->
Mmin
,
p
))
/
p
;
}
for
(
i
=
0
;
i
<
n
-
1
;
i
++
)
{
cte
=
cte
*
pow
(
Cp
->
ms
[
i
],(
Cp
->
as
[
i
]
-
Cp
->
as
[
i
+
1
]
));
p
=
Cp
->
as
[
i
+
1
]
+
1
;
integral
=
integral
+
cte
*
(
pow
(
Cp
->
ms
[
i
+
1
],
p
)
-
pow
(
Cp
->
ms
[
i
],
p
))
/
p
;
}
if
(
Cp
->
Mmax
>
Cp
->
ms
[
-
1
])
{
cte
=
cte
*
pow
(
Cp
->
ms
[
n
-
1
]
,
(
Cp
->
as
[
n
-
1
]
-
Cp
->
as
[
n
]
)
);
p
=
Cp
->
as
[
n
]
+
1
;
integral
=
integral
+
cte
*
(
pow
(
Cp
->
Mmax
,
p
)
-
pow
(
Cp
->
ms
[
n
-
1
],
p
))
/
p
;
}
/* compute all b */
Cp
->
bs
[
0
]
=
1.
/
integral
;
for
(
i
=
0
;
i
<
n
;
i
++
)
{
Cp
->
bs
[
i
+
1
]
=
Cp
->
bs
[
i
]
*
pow
(
Cp
->
ms
[
i
],(
Cp
->
as
[
i
]
-
Cp
->
as
[
i
+
1
]
));
}
}
if
(
verbose
&&
ThisTask
==
0
)
{
printf
(
"-- bs --
\n
"
);
for
(
i
=
0
;
i
<
n
+
1
;
i
++
)
printf
(
"%g "
,
Cp
->
bs
[
i
]);
printf
(
"
\n
"
);
}
mmin
=
Cp
->
Mmin
/
All
.
UnitMass_in_g
*
SOLAR_MASS
;
/* in mass unit */
mmax
=
Cp
->
Mmax
/
All
.
UnitMass_in_g
*
SOLAR_MASS
;
/* in mass unit */
Cp
->
imf_Ntot
=
get_imf_N
(
mmin
,
mmax
)
*
SOLAR_MASS
/
All
.
UnitMass_in_g
;
/* init fs : mass fraction at ms */
if
(
n
>
0
)
{
for
(
i
=
0
;
i
<
n
+
1
;
i
++
)
{
mmax
=
Cp
->
ms
[
i
]
/
All
.
UnitMass_in_g
*
SOLAR_MASS
;
/* in mass unit */
Cp
->
fs
[
i
]
=
SOLAR_MASS
/
All
.
UnitMass_in_g
*
get_imf_N
(
mmin
,
mmax
)
/
Cp
->
imf_Ntot
;
}
}
}
/*
This function init the chime parameters
*/
void
init_chimie
(
void
)
{
int
i
,
nf
;
double
u_lt
;
double
UnitLength_in_kpc
;
double
UnitMass_in_Msol
;
char
filename
[
500
];
char
ext
[
100
];
/* check some flags */
#ifndef COSMICTIME
if
(
All
.
ComovingIntegrationOn
)
{
if
(
ThisTask
==
0
)
printf
(
"Code wasn't compiled with COSMICTIME support enabled!
\n
"
);
endrun
(
-
88800
);
}
#endif
UnitLength_in_kpc
=
All
.
UnitLength_in_cm
/
KPC_IN_CM
;
UnitMass_in_Msol
=
All
.
UnitMass_in_g
/
SOLAR_MASS
;
//u_lt = -log10( 4.7287e11*sqrt(pow(UnitLength_in_kpc,3)/UnitMass_in_Msol));
/*Sat Dec 25 23:27:10 CET 2010 */
u_lt
=
-
log10
(
All
.
UnitTime_in_Megayears
*
1e6
);
allocate_chimie
();
for
(
nf
=
0
;
nf
<
All
.
ChimieNumberOfParameterFiles
;
nf
++
)
{
if
(
All
.
ChimieNumberOfParameterFiles
==
1
)
sprintf
(
filename
,
"%s"
,
All
.
ChimieParameterFile
);
else
sprintf
(
filename
,
"%s.%d"
,
All
.
ChimieParameterFile
,
nf
);
read_chimie
(
filename
,
nf
);
/* set the table */
set_table
(
nf
);
/* Conversion into program time unit */
Cp
->
coeff_z
[
2
][
2
]
=
Cp
->
coeff_z
[
2
][
2
]
+
u_lt
;
for
(
i
=
0
;
i
<
3
;
i
++
)
Cp
->
coeff_z
[
1
][
i
]
=
Cp
->
coeff_z
[
1
][
i
]
/
2.0
;
/* init imf parameters */
init_imf
();
/* init SNII parameters */
if
(
Cp
->
n
==
0
)
{
//Cp->SNII_cte[0] = Cp->bs[0]/Cp->as[0];
Cp
->
SNII_cte
=
Cp
->
bs
[
0
]
/
Cp
->
as
[
0
];
Cp
->
SNII_a
=
Cp
->
as
[
0
];
}
else
{
//for (i=0;i<Cp->n+1;i++) /* if multiple power law in the SNII mass range */
// Cp->SNII_cte[i] = Cp->bs[i]/Cp->as[i];
Cp
->
SNII_cte
=
Cp
->
bs
[
Cp
->
n
]
/
Cp
->
as
[
Cp
->
n
];
Cp
->
SNII_a
=
Cp
->
as
[
Cp
->
n
];
}
/* init SNIa parameters */
Cp
->
SNIa_a1
=
Cp
->
SNIa_a
;
Cp
->
SNIa_b1
=
(
Cp
->
SNIa_a1
+
1
)
/
(
pow
(
Cp
->
SNIa_Mdu1
,
Cp
->
SNIa_a1
+
1
)
-
pow
(
Cp
->
SNIa_Mdl1
,
Cp
->
SNIa_a1
+
1
));
Cp
->
SNIa_cte1
=
Cp
->
SNIa_b1
/
Cp
->
SNIa_a1
;
Cp
->
SNIa_a2
=
Cp
->
SNIa_a
;
Cp
->
SNIa_b2
=
(
Cp
->
SNIa_a2
+
1
)
/
(
pow
(
Cp
->
SNIa_Mdu2
,
Cp
->
SNIa_a2
+
1
)
-
pow
(
Cp
->
SNIa_Mdl2
,
Cp
->
SNIa_a2
+
1
));
Cp
->
SNIa_cte2
=
Cp
->
SNIa_b2
/
Cp
->
SNIa_a2
;
/* init SNII parameters */
if
(
Cp
->
n
==
0
)
{
Cp
->
SNIa_cte
=
Cp
->
bs
[
0
]
/
Cp
->
as
[
0
];
Cp
->
SNIa_a
=
Cp
->
as
[
0
];
}
else
{
Cp
->
SNIa_cte
=
Cp
->
bs
[
Cp
->
n
]
/
Cp
->
as
[
Cp
->
n
];
Cp
->
SNIa_a
=
Cp
->
as
[
Cp
->
n
];
}
Cp
->
SNII_Mmax
=
Cp
->
Mmax
;
for
(
i
=
0
;
i
<
Cp
->
nelts
+
2
;
i
++
)
Elt
[
i
].
Mmin
=
log10
(
Elt
[
i
].
Mmin
);
/* output info */
if
(
verbose
&&
ThisTask
==
0
)
{
printf
(
"-- SNII_cte --
\n
"
);
//for (i=0;i<Cp->n+1;i++)
// printf("%g ",Cp->SNII_cte[i]);
printf
(
"%g "
,
Cp
->
SNII_cte
);
printf
(
"
\n
"
);
}
/* check that the masses are higher than the last IMF elbow */
if
(
Cp
->
n
>
0
)
{
if
(
Cp
->
SNIa_Mpl
<
Cp
->
ms
[
Cp
->
n
-
1
])
{
printf
(
"
\n
SNIa_Mpl = %g < ms[n-1] = %g !!!
\n\n
"
,
Cp
->
SNIa_Mpl
,
Cp
->
ms
[
Cp
->
n
-
1
]);
printf
(
"This is not supported by the current implementation !!!
\n
"
);
endrun
(
88801
);
}
if
(
Cp
->
SNIa_Mpu
<
Cp
->
ms
[
Cp
->
n
-
1
])
{
printf
(
"
\n
SNIa_Mpu = %g < ms[n-1] = %g !!!
\n\n
"
,
Cp
->
SNIa_Mpu
,
Cp
->
ms
[
Cp
->
n
-
1
]);
printf
(
"This is not supported by the current implementation !!!
\n
"
);
endrun
(
88802
);
}
if
(
Cp
->
SNII_Mmin
<
Cp
->
ms
[
Cp
->
n
-
1
])
{
printf
(
"
\n
SNII_Mmin = %g < ms[n-1] = %g !!!
\n\n
"
,
Cp
->
SNII_Mmin
,
Cp
->
ms
[
Cp
->
n
-
1
]);
printf
(
"This is not supported by the current implementation !!!
\n
"
);
endrun
(
88803
);
}
if
(
Cp
->
SNII_Mmax
<
Cp
->
ms
[
Cp
->
n
-
1
])
{
printf
(
"
\n
SNII_Mmax = %g < ms[n-1] = %g !!!
\n\n
"
,
Cp
->
SNII_Mmax
,
Cp
->
ms
[
Cp
->
n
-
1
]);
printf
(
"This is not supported by the current implementation !!!
\n
"
);
endrun
(
88804
);
}
}
}
}
void
check_chimie
(
void
)
{
int
i
;
printf
(
"(Taks=%d) Number of elts : %d
\n
"
,
ThisTask
,
Cp
->
nelts
);
for
(
i
=
2
;
i
<
Cp
->
nelts
+
2
;
i
++
)
printf
(
"%s "
,
&
Elt
[
i
].
label
);
printf
(
"
\n
"
);
/* check number of elements */
if
(
NELEMENTS
!=
Cp
->
nelts
)
{
printf
(
"(Taks=%d) NELEMENTS (=%d) != Cp->nelts (=%d) : please check !!!
\n\n
"
,
ThisTask
,
NELEMENTS
,
Cp
->
nelts
);
endrun
(
88807
);
}
/* check that iron is the first element */
if
((
strcmp
(
"Fe"
,
Elt
[
2
].
label
))
!=
0
)
{
printf
(
"(Taks=%d) first element (=%s) is not %s !!!
\n\n
"
,
ThisTask
,
Elt
[
2
].
label
,
FIRST_ELEMENT
);
endrun
(
88808
);
}
}
int
get_nelts
()
{
return
Cp
->
nelts
;
}
float
get_SolarAbundance
(
i
)
{
return
Elt
[
i
+
2
].
SolarAbundance
;
}
char
*
get_Element
(
i
)
{
return
Elt
[
i
+
2
].
label
;
}
double
star_lifetime
(
double
z
,
double
m
)
{
/* z is the mass fraction of metals, ie, the metallicity */
/* m is the stellar mass in code unit */
/* Return t in code time unit */
int
i
;
double
a
,
b
,
c
;
double
coeff
[
3
];
double
logm
,
twologm
,
logm2
,
time
;
/* convert m in msol */
m
=
m
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
for
(
i
=
0
;
i
<
3
;
i
++
)
coeff
[
i
]
=
(
Cp
->
coeff_z
[
i
][
0
]
*
z
+
Cp
->
coeff_z
[
i
][
1
]
)
*
z
+
Cp
->
coeff_z
[
i
][
2
];
a
=
coeff
[
0
];
b
=
coeff
[
1
];
c
=
coeff
[
2
];
logm
=
log10
(
m
);
twologm
=
2.0
*
logm
;
logm2
=
logm
*
logm
;
time
=
pow
(
10.
,(
a
*
logm2
+
b
*
twologm
+
c
));
return
time
;
}
double
star_mass_from_age
(
double
z
,
double
t
)
{
/* z is the mass fraction of metals, ie, the metallicity */
/* t is the star life time */
/* return the stellar mass (in code unit) that has a lifetime equal to t */
/* this is the inverse of star_lifetime */
int
i
;
double
a
,
b
,
c
;
double
coeff
[
3
];
double
m
;
for
(
i
=
0
;
i
<
3
;
i
++
)
coeff
[
i
]
=
(
Cp
->
coeff_z
[
i
][
0
]
*
z
+
Cp
->
coeff_z
[
i
][
1
]
)
*
z
+
Cp
->
coeff_z
[
i
][
2
];
a
=
coeff
[
0
];
b
=
coeff
[
1
];
c
=
coeff
[
2
];
m
=
-
(
b
+
sqrt
(
b
*
b
-
a
*
(
c
-
log10
(
t
))))
/
a
;
m
=
pow
(
10
,
m
);
/* here, m is in solar mass */
m
=
m
*
SOLAR_MASS
/
All
.
UnitMass_in_g
;
/* Msol to mass unit */
return
m
;
}
/****************************************************************************************/
/*
/* Supernova rate : number of supernova per mass unit
/*
/****************************************************************************************/
double
SNII_rate
(
double
m1
,
double
m2
)
{
/*
masses in code unit
*/
double
RSNII
;
double
md
,
mu
;
RSNII
=
0.0
;
/* convert m in msol */
m1
=
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
m2
=
m2
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
/* (1) find md, mu */
md
=
dmax
(
m1
,
Cp
->
SNII_Mmin
);
mu
=
dmin
(
m2
,
Cp
->
SNII_Mmax
);
if
(
mu
<=
md
)
/* no SNII in that mass range */
return
0.0
;
RSNII
=
Cp
->
SNII_cte
*
(
pow
(
mu
,
Cp
->
SNII_a
)
-
pow
(
md
,
Cp
->
SNII_a
));
/* number per solar mass */
/* convert in number per solar mass to number per mass unit */
RSNII
=
RSNII
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
return
RSNII
;
}
double
SNIa_rate
(
double
m1
,
double
m2
)
{
/*
masses in code unit
*/
double
RSNIa
;
double
md
,
mu
;
RSNIa
=
0.0
;
/* convert m in msol */
m1
=
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
m2
=
m2
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
/* RG contribution */
md
=
dmax
(
m1
,
Cp
->
SNIa_Mdl1
);
mu
=
dmin
(
m2
,
Cp
->
SNIa_Mdu1
);
if
(
md
<
mu
)
RSNIa
=
RSNIa
+
Cp
->
SNIa_bb1
*
Cp
->
SNIa_cte1
*
(
pow
(
mu
,
Cp
->
SNIa_a1
)
-
pow
(
md
,
Cp
->
SNIa_a1
));
/* MS contribution */
md
=
dmax
(
m1
,
Cp
->
SNIa_Mdl2
);
mu
=
dmin
(
m2
,
Cp
->
SNIa_Mdu2
);
if
(
md
<
mu
)
RSNIa
=
RSNIa
+
Cp
->
SNIa_bb2
*
Cp
->
SNIa_cte2
*
(
pow
(
mu
,
Cp
->
SNIa_a2
)
-
pow
(
md
,
Cp
->
SNIa_a2
));
/* WD contribution */
md
=
dmax
(
m1
,
Cp
->
SNIa_Mpl
);
/* select stars that have finished their life -> WD */
mu
=
Cp
->
SNIa_Mpu
;
/* no upper bond */
if
(
mu
<=
md
)
/* no SNIa in that mass range */
return
0.0
;
RSNIa
=
RSNIa
*
Cp
->
SNIa_cte
*
(
pow
(
mu
,
Cp
->
SNIa_a
)
-
pow
(
md
,
Cp
->
SNIa_a
));
/* number per solar mass */
/* convert in number per solar mass to number per mass unit */
RSNIa
=
RSNIa
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
return
RSNIa
;
}
void
SNII_mass_ejection
(
double
m1
,
double
m2
)
{
double
l1
,
l2
;
int
i1
,
i2
,
i1p
,
i2p
,
j
;
double
f1
,
f2
;
double
v1
,
v2
;
/* convert m in msol */
m1
=
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
m2
=
m2
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
j
=
0
;
l1
=
(
log10
(
m1
)
-
Elt
[
j
].
Mmin
)
/
Elt
[
j
].
Step
;
l2
=
(
log10
(
m2
)
-
Elt
[
j
].
Mmin
)
/
Elt
[
j
].
Step
;
if
(
l1
<
0.0
)
l1
=
0.0
;
if
(
l2
<
0.0
)
l2
=
0.0
;
i1
=
(
int
)
l1
;
i2
=
(
int
)
l2
;
i1p
=
i1
+
1
;
i2p
=
i2
+
1
;
f1
=
l1
-
i1
;
f2
=
l2
-
i2
;
/* check (yr) */
if
(
i1
<
0
)
i1
=
0
;
if
(
i2
<
0
)
i2
=
0
;
/* --------- TOTAL GAS ---------- */
j
=
0
;
v1
=
f1
*
(
Elt
[
j
].
Array
[
i1p
]
-
Elt
[
j
].
Array
[
i1
]
)
+
Elt
[
j
].
Array
[
i1
];
v2
=
f2
*
(
Elt
[
j
].
Array
[
i2p
]
-
Elt
[
j
].
Array
[
i2
]
)
+
Elt
[
j
].
Array
[
i2
];
MassFracSNII
[
j
]
=
v2
-
v1
;
/* --------- He core therm ---------- */
j
=
1
;
v1
=
f1
*
(
Elt
[
j
].
Array
[
i1p
]
-
Elt
[
j
].
Array
[
i1
]
)
+
Elt
[
j
].
Array
[
i1
];
v2
=
f2
*
(
Elt
[
j
].
Array
[
i2p
]
-
Elt
[
j
].
Array
[
i2
]
)
+
Elt
[
j
].
Array
[
i2
];
MassFracSNII
[
j
]
=
v2
-
v1
;
/* ---------------------------- */
/* --------- Metals ---------- */
/* ---------------------------- */
j
=
2
;
l1
=
(
log10
(
m1
)
-
Elt
[
j
].
Mmin
)
/
Elt
[
j
].
Step
;
l2
=
(
log10
(
m2
)
-
Elt
[
j
].
Mmin
)
/
Elt
[
j
].
Step
;
if
(
l1
<
0.0
)
l1
=
0.0
;
if
(
l2
<
0.0
)
l2
=
0.0
;
i1
=
(
int
)
l1
;
i2
=
(
int
)
l2
;
i1p
=
i1
+
1
;
i2p
=
i2
+
1
;
f1
=
l1
-
i1
;
f2
=
l2
-
i2
;
/* check (yr) */
if
(
i1
<
0
)
i1
=
0
;
if
(
i2
<
0
)
i2
=
0
;
for
(
j
=
2
;
j
<
Cp
->
nelts
+
2
;
j
++
)
{
v1
=
f1
*
(
Elt
[
j
].
Array
[
i1p
]
-
Elt
[
j
].
Array
[
i1
]
)
+
Elt
[
j
].
Array
[
i1
];
v2
=
f2
*
(
Elt
[
j
].
Array
[
i2p
]
-
Elt
[
j
].
Array
[
i2
]
)
+
Elt
[
j
].
Array
[
i2
];
MassFracSNII
[
j
]
=
v2
-
v1
;
}
}
void
SNII_single_mass_ejection
(
double
m1
)
{
double
l1
;
int
i1
,
i1p
,
j
;
double
f1
;
double
v1
;
/* convert m in msol */
m1
=
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
;
j
=
0
;
l1
=
(
log10
(
m1
)
-
Elt
[
j
].
Mmin
)
/
Elt
[
j
].
Step
;
if
(
l1
<
0.0
)
l1
=
0.0
;
i1
=
(
int
)
l1
;
i1p
=
i1
+
1
;
f1
=
l1
-
i1
;
/* check (yr) */
if
(
i1
<
0
)
i1
=
0
;
/* --------- TOTAL GAS ---------- */
j
=
0
;
v1
=
f1
*
(
Elt
[
j
].
Metal
[
i1p
]
-
Elt
[
j
].
Metal
[
i1
]
)
+
Elt
[
j
].
Metal
[
i1
];
SingleMassFracSNII
[
j
]
=
v1
;
/* --------- He core therm ---------- */
j
=
1
;
v1
=
f1
*
(
Elt
[
j
].
Metal
[
i1p
]
-
Elt
[
j
].
Metal
[
i1
]
)
+
Elt
[
j
].
Metal
[
i1
];
SingleMassFracSNII
[
j
]
=
v1
;
/* ---------------------------- */
/* --------- Metals ---------- */
/* ---------------------------- */
j
=
2
;
l1
=
(
log10
(
m1
)
-
Elt
[
j
].
Mmin
)
/
Elt
[
j
].
Step
;
if
(
l1
<
0.0
)
l1
=
0.0
;
i1
=
(
int
)
l1
;
i1p
=
i1
+
1
;
f1
=
l1
-
i1
;
/* check (yr) */
if
(
i1
<
0
)
i1
=
0
;
for
(
j
=
2
;
j
<
Cp
->
nelts
+
2
;
j
++
)
{
v1
=
f1
*
(
Elt
[
j
].
Metal
[
i1p
]
-
Elt
[
j
].
Metal
[
i1
]
)
+
Elt
[
j
].
Metal
[
i1
];
SingleMassFracSNII
[
j
]
=
v1
;
}
}
void
Total_mass_ejection
(
double
m1
,
double
m2
,
double
M0
,
double
*
z
)
{
int
j
;
double
NSNIa
;
/* compute SNII mass ejection -> MassFracSNII */
SNII_mass_ejection
(
m1
,
m2
);
/* number of SNIa per mass unit between time and time+dt */
NSNIa
=
SNIa_rate
(
m1
,
m2
)
*
M0
;
/* number of SNII per mass unit between time and time+dt */
//NSNII = SNII_rate(m1,m2)*M0; /* useless (only for energy) */
/* total ejected gas mass */
EjectedMass
[
0
]
=
M0
*
MassFracSNII
[
0
]
+
Cp
->
Mco
/
All
.
UnitMass_in_g
*
SOLAR_MASS
*
NSNIa
;
/* ejected mass per element */
for
(
j
=
2
;
j
<
Cp
->
nelts
+
2
;
j
++
)
EjectedMass
[
j
]
=
M0
*
(
MassFracSNII
[
j
]
+
z
[
j
-
2
]
*
MassFracSNII
[
1
])
+
NSNIa
*
Elt
[
j
].
MSNIa
/
All
.
UnitMass_in_g
*
SOLAR_MASS
;
/* not used */
EjectedMass
[
1
]
=
-
1
;
}
void
Total_single_mass_ejection
(
double
m1
,
double
*
z
)
{
/*
!!! we do not take into account SNIa
*/
int
j
;
float
M0
;
M0
=
m1
;
/* compute SNII mass ejection -> SingleMassFracSNII */
SNII_single_mass_ejection
(
m1
);
/* total ejected gas mass */
SingleEjectedMass
[
0
]
=
M0
*
SingleMassFracSNII
[
0
];
/* + Cp->Mco/All.UnitMass_in_g*SOLAR_MASS * NSNIa; */
/* ejected mass per element */
for
(
j
=
2
;
j
<
Cp
->
nelts
+
2
;
j
++
)
SingleEjectedMass
[
j
]
=
M0
*
(
SingleMassFracSNII
[
j
]
+
z
[
j
-
2
]
*
SingleMassFracSNII
[
1
]);
/* + NSNIa* Elt[j].MSNIa/All.UnitMass_in_g*SOLAR_MASS; */
/* not used */
SingleEjectedMass
[
1
]
=
-
1
;
}
/**********************************************************************************************
END OF CHIMIE FUNCTIONS
**********************************************************************************************/
#if defined(CHIMIE_THERMAL_FEEDBACK) && defined(CHIMIE_COMPUTE_THERMAL_FEEDBACK_ENERGY)
void
chimie_compute_energy_int
(
int
mode
)
{
int
i
;
double
DeltaEgyInt
;
double
Tot_DeltaEgyInt
;
DeltaEgyInt
=
0
;
Tot_DeltaEgyInt
=
0
;
if
(
mode
==
1
)
{
LocalSysState
.
EnergyInt1
=
0
;
LocalSysState
.
EnergyInt2
=
0
;
}
for
(
i
=
0
;
i
<
N_gas
;
i
++
)
{
if
(
P
[
i
].
Type
==
0
)
{
if
(
mode
==
1
)
LocalSysState
.
EnergyInt1
+=
P
[
i
].
Mass
*
SphP
[
i
].
EntropyPred
/
(
GAMMA_MINUS1
)
*
pow
(
SphP
[
i
].
Density
*
a3inv
,
GAMMA_MINUS1
);
else
LocalSysState
.
EnergyInt2
+=
P
[
i
].
Mass
*
SphP
[
i
].
EntropyPred
/
(
GAMMA_MINUS1
)
*
pow
(
SphP
[
i
].
Density
*
a3inv
,
GAMMA_MINUS1
);
}
}
if
(
mode
==
2
)
{
DeltaEgyInt
=
LocalSysState
.
EnergyInt2
-
LocalSysState
.
EnergyInt1
;
MPI_Reduce
(
&
DeltaEgyInt
,
&
Tot_DeltaEgyInt
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
LocalSysState
.
EnergyThermalFeedback
-=
DeltaEgyInt
;
}
}
#endif
#if defined(CHIMIE_KINETIC_FEEDBACK) && defined(CHIMIE_COMPUTE_KINETIC_FEEDBACK_ENERGY)
void
chimie_compute_energy_kin
(
int
mode
)
{
int
i
;
double
DeltaEgyKin
;
double
Tot_DeltaEgyKin
;
DeltaEgyKin
=
0
;
Tot_DeltaEgyKin
=
0
;
if
(
mode
==
1
)
{
LocalSysState
.
EnergyKin1
=
0
;
LocalSysState
.
EnergyKin2
=
0
;
}
for
(
i
=
0
;
i
<
N_gas
;
i
++
)
{
if
(
P
[
i
].
Type
==
0
)
{
if
(
mode
==
1
)
LocalSysState
.
EnergyKin1
+=
0.5
*
P
[
i
].
Mass
*
(
P
[
i
].
Vel
[
0
]
*
P
[
i
].
Vel
[
0
]
+
P
[
i
].
Vel
[
1
]
*
P
[
i
].
Vel
[
1
]
+
P
[
i
].
Vel
[
2
]
*
P
[
i
].
Vel
[
2
]);
else
LocalSysState
.
EnergyKin2
+=
0.5
*
P
[
i
].
Mass
*
(
P
[
i
].
Vel
[
0
]
*
P
[
i
].
Vel
[
0
]
+
P
[
i
].
Vel
[
1
]
*
P
[
i
].
Vel
[
1
]
+
P
[
i
].
Vel
[
2
]
*
P
[
i
].
Vel
[
2
]);
}
}
if
(
mode
==
2
)
{
DeltaEgyKin
=
LocalSysState
.
EnergyKin2
-
LocalSysState
.
EnergyKin1
;
MPI_Reduce
(
&
DeltaEgyKin
,
&
Tot_DeltaEgyKin
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
LocalSysState
.
EnergyKineticFeedback
-=
DeltaEgyKin
;
}
}
#endif
#ifdef CHIMIE_THERMAL_FEEDBACK
void
chimie_apply_thermal_feedback
(
void
)
{
int
i
;
double
EgySpec
,
NewEgySpec
,
DeltaEntropy
;
for
(
i
=
0
;
i
<
N_gas
;
i
++
)
{
if
(
P
[
i
].
Type
==
0
)
{
if
(
SphP
[
i
].
DeltaEgySpec
>
0
)
{
/* spec energy at current step */
EgySpec
=
SphP
[
i
].
EntropyPred
/
GAMMA_MINUS1
*
pow
(
SphP
[
i
].
Density
*
a3inv
,
GAMMA_MINUS1
);
/* new egyspec */
NewEgySpec
=
EgySpec
+
SphP
[
i
].
DeltaEgySpec
;
LocalSysState
.
EnergyThermalFeedback
-=
SphP
[
i
].
DeltaEgySpec
*
P
[
i
].
Mass
;
/* new entropy */
DeltaEntropy
=
GAMMA_MINUS1
*
NewEgySpec
/
pow
(
SphP
[
i
].
Density
*
a3inv
,
GAMMA_MINUS1
)
-
SphP
[
i
].
EntropyPred
;
SphP
[
i
].
EntropyPred
+=
DeltaEntropy
;
SphP
[
i
].
Entropy
+=
DeltaEntropy
;
/* reset variable */
SphP
[
i
].
DeltaEgySpec
=
0
;
/* recode time */
SphP
[
i
].
ThermalTime
=
All
.
Time
;
}
}
}
}
#endif
#ifdef CHIMIE_KINETIC_FEEDBACK
void
chimie_apply_wind
(
void
)
{
/* apply wind */
int
i
;
double
e1
,
e2
;
double
phi
,
costh
,
sinth
,
vx
,
vy
,
vz
;
for
(
i
=
0
;
i
<
N_gas
;
i
++
)
{
if
(
P
[
i
].
Type
==
0
)
{
if
(
SphP
[
i
].
WindFlag
)
{
phi
=
get_ChimieKineticFeedback_random_number
(
P
[
i
].
ID
)
*
PI
*
2.
;
costh
=
1.
-
2.
*
get_ChimieKineticFeedback_random_number
(
P
[
i
].
ID
+
1
);
sinth
=
sqrt
(
1.
-
pow
(
costh
,
2
));
vx
=
All
.
ChimieWindSpeed
*
sinth
*
cos
(
phi
);
vy
=
All
.
ChimieWindSpeed
*
sinth
*
sin
(
phi
);
vz
=
All
.
ChimieWindSpeed
*
costh
;
e1
=
0.5
*
P
[
i
].
Mass
*
(
SphP
[
i
].
VelPred
[
0
]
*
SphP
[
i
].
VelPred
[
0
]
+
SphP
[
i
].
VelPred
[
1
]
*
SphP
[
i
].
VelPred
[
1
]
+
SphP
[
i
].
VelPred
[
2
]
*
SphP
[
i
].
VelPred
[
2
]);
P
[
i
].
Vel
[
0
]
+=
vx
;
P
[
i
].
Vel
[
1
]
+=
vy
;
P
[
i
].
Vel
[
2
]
+=
vz
;
SphP
[
i
].
VelPred
[
0
]
+=
vx
;
SphP
[
i
].
VelPred
[
1
]
+=
vy
;
SphP
[
i
].
VelPred
[
2
]
+=
vz
;
e2
=
0.5
*
P
[
i
].
Mass
*
(
SphP
[
i
].
VelPred
[
0
]
*
SphP
[
i
].
VelPred
[
0
]
+
SphP
[
i
].
VelPred
[
1
]
*
SphP
[
i
].
VelPred
[
1
]
+
SphP
[
i
].
VelPred
[
2
]
*
SphP
[
i
].
VelPred
[
2
]);
LocalSysState
.
EnergyKineticFeedback
-=
e2
-
e1
;
SphP
[
i
].
WindFlag
=
0
;
}
}
}
}
#endif
/*! This function is the driver routine for the calculation of chemical evolution
*/
void
chimie
(
void
)
{
double
t0
,
t1
;
t0
=
second
();
/* measure the time for the full chimie computation */
if
(
ThisTask
==
0
)
printf
(
"Start Chimie computation.
\n
"
);
/* apply thermal feedback on selected particles */
#ifdef CHIMIE_THERMAL_FEEDBACK
chimie_apply_thermal_feedback
();
#endif
/* apply wind on selected particles */
#ifdef CHIMIE_KINETIC_FEEDBACK
chimie_apply_wind
();
#endif
stars_density
();
/* compute density */
do_chimie
();
/* chimie */
if
(
ThisTask
==
0
)
printf
(
"Chimie computation done.
\n
"
);
t1
=
second
();
All
.
CPU_Chimie
+=
timediff
(
t0
,
t1
);
}
/*! This function is the driver routine for the calculation of chemical evolution
*/
void
do_chimie
(
void
)
{
long
long
ntot
,
ntotleft
;
int
i
,
j
,
k
,
n
,
m
,
ngrp
,
maxfill
,
source
,
ndone
;
int
*
nbuffer
,
*
noffset
,
*
nsend_local
,
*
nsend
,
*
numlist
,
*
ndonelist
;
int
level
,
sendTask
,
recvTask
,
nexport
,
place
;
double
tstart
,
tend
,
sumt
,
sumcomm
;
double
timecomp
=
0
,
timecommsumm
=
0
,
timeimbalance
=
0
,
sumimbalance
;
int
flag_chimie
;
MPI_Status
status
;
int
do_it
;
int
Ti0
,
Ti1
,
Ti2
;
double
t1
,
t2
,
t01
,
t02
;
double
tmin
,
tmax
;
double
minlivetime
,
maxlivetime
;
double
m1
,
m2
,
M0
;
double
NSNIa
,
NSNII
;
double
NSNIa_tot
,
NSNII_tot
,
NSNIa_totlocal
,
NSNII_totlocal
;
double
EgySN
,
EgySNlocal
;
double
EgySNThermal
,
EgySNKinetic
;
int
Nchim
,
Nchimlocal
;
int
Nwind
,
Nwindlocal
;
int
Nflag
,
Nflaglocal
;
int
Noldwind
,
Noldwindlocal
;
double
metals
[
NELEMENTS
];
double
FeH
;
float
MinRelMass
=
1e-3
;
#ifdef PERIODIC
boxSize
=
All
.
BoxSize
;
boxHalf
=
0.5
*
All
.
BoxSize
;
#ifdef LONG_X
boxHalf_X
=
boxHalf
*
LONG_X
;
boxSize_X
=
boxSize
*
LONG_X
;
#endif
#ifdef LONG_Y
boxHalf_Y
=
boxHalf
*
LONG_Y
;
boxSize_Y
=
boxSize
*
LONG_Y
;
#endif
#ifdef LONG_Z
boxHalf_Z
=
boxHalf
*
LONG_Z
;
boxSize_Z
=
boxSize
*
LONG_Z
;
#endif
#endif
#ifdef COMPUTE_VELOCITY_DISPERSION
double
v1m
,
v2m
;
#endif
if
(
All
.
ComovingIntegrationOn
)
{
/* Factors for comoving integration of hydro */
hubble_a
=
All
.
Omega0
/
(
All
.
Time
*
All
.
Time
*
All
.
Time
)
+
(
1
-
All
.
Omega0
-
All
.
OmegaLambda
)
/
(
All
.
Time
*
All
.
Time
)
+
All
.
OmegaLambda
;
hubble_a
=
All
.
Hubble
*
sqrt
(
hubble_a
);
hubble_a2
=
All
.
Time
*
All
.
Time
*
hubble_a
;
fac_mu
=
pow
(
All
.
Time
,
3
*
(
GAMMA
-
1
)
/
2
)
/
All
.
Time
;
fac_egy
=
pow
(
All
.
Time
,
3
*
(
GAMMA
-
1
));
fac_vsic_fix
=
hubble_a
*
pow
(
All
.
Time
,
3
*
GAMMA_MINUS1
);
a3inv
=
1
/
(
All
.
Time
*
All
.
Time
*
All
.
Time
);
atime
=
All
.
Time
;
#ifdef FEEDBACK
fac_pow
=
fac_egy
*
atime
*
atime
;
#endif
}
else
{
hubble_a
=
hubble_a2
=
atime
=
fac_mu
=
fac_vsic_fix
=
a3inv
=
fac_egy
=
1.0
;
#ifdef FEEDBACK
fac_pow
=
1.0
;
#endif
}
/* `NumStUpdate' gives the number of particles on this processor that want a chimie computation */
for
(
n
=
0
,
NumStUpdate
=
0
;
n
<
N_gas
+
N_stars
;
n
++
)
{
if
(
P
[
n
].
Ti_endstep
==
All
.
Ti_Current
)
if
(
P
[
n
].
Type
==
ST
)
{
m
=
P
[
n
].
StPIdx
;
if
(
(
P
[
n
].
Mass
/
StP
[
m
].
InitialMass
)
>
MinRelMass
)
NumStUpdate
++
;
}
if
(
P
[
n
].
Type
==
0
)
SphP
[
n
].
dMass
=
0.
;
}
numlist
=
malloc
(
NTask
*
sizeof
(
int
)
*
NTask
);
MPI_Allgather
(
&
NumStUpdate
,
1
,
MPI_INT
,
numlist
,
1
,
MPI_INT
,
MPI_COMM_WORLD
);
for
(
i
=
0
,
ntot
=
0
;
i
<
NTask
;
i
++
)
ntot
+=
numlist
[
i
];
free
(
numlist
);
noffset
=
malloc
(
sizeof
(
int
)
*
NTask
);
/* offsets of bunches in common list */
nbuffer
=
malloc
(
sizeof
(
int
)
*
NTask
);
nsend_local
=
malloc
(
sizeof
(
int
)
*
NTask
);
nsend
=
malloc
(
sizeof
(
int
)
*
NTask
*
NTask
);
ndonelist
=
malloc
(
sizeof
(
int
)
*
NTask
);
i
=
0
;
/* first gas particle, because stars may be hidden among gas particles */
ntotleft
=
ntot
;
/* particles left for all tasks together */
NSNIa_tot
=
0
;
NSNII_tot
=
0
;
NSNIa_totlocal
=
0
;
NSNII_totlocal
=
0
;
EgySN
=
0
;
EgySNlocal
=
0
;
Nchimlocal
=
0
;
Nchim
=
0
;
Nwindlocal
=
0
;
Nwind
=
0
;
Noldwindlocal
=
0
;
Noldwind
=
0
;
Nflaglocal
=
0
;
Nflag
=
0
;
while
(
ntotleft
>
0
)
{
for
(
j
=
0
;
j
<
NTask
;
j
++
)
nsend_local
[
j
]
=
0
;
/* do local particles and prepare export list */
tstart
=
second
();
for
(
nexport
=
0
,
ndone
=
0
;
i
<
N_gas
+
N_stars
&&
nexport
<
All
.
BunchSizeChimie
-
NTask
;
i
++
)
{
/* only active particles and stars */
if
((
P
[
i
].
Ti_endstep
==
All
.
Ti_Current
)
&&
(
P
[
i
].
Type
==
ST
))
{
if
(
P
[
i
].
Type
!=
ST
)
{
printf
(
"P[i].Type != ST, we better stop.
\n
"
);
printf
(
"N_gas=%d (type=%d) i=%d (type=%d)
\n
"
,
N_gas
,
P
[
N_gas
].
Type
,
i
,
P
[
i
].
Type
);
printf
(
"Please, check that you do not use PEANOHILBERT
\n
"
);
endrun
(
777001
);
}
m
=
P
[
i
].
StPIdx
;
if
(
(
P
[
i
].
Mass
/
StP
[
m
].
InitialMass
)
>
MinRelMass
)
{
flag_chimie
=
0
;
/******************************************/
/* do chimie */
/******************************************/
/*****************************************************/
/* look if a SN may have explode during the last step
/*****************************************************/
/***********************************************/
/***********************************************/
/* set the right table base of the metallicity */
set_table
(
0
);
//FeH = log10( (StP[m].Metal[FE]/get_SolarAbundance(FE)) + 1.e-20 );
//if (FeH<-3)
// set_table(1);
//else
// set_table(0);
//if (P[i].ID==65546)
// {
// printf("(%d) %g the particle 65546 FeH=%g metalFe=%g Mmin=%g Mmax=%g n=%d\n",ThisTask,All.Time,FeH,StP[m].Metal[FE],Cp->Mmin,Cp->Mmax,Cp->n);
// }
/*
Cp->Mmin
Cp->Mmax
Cp->n
Cp->ms[]
Cp->as[]
Cp->SNIa_cte
Cp->SNIa_a
Cp->SNIa_Mdl1
Cp->SNIa_Mdu1
Cp->SNIa_bb1
Cp->SNIa_cte1
Cp->SNIa_a1
Cp->SNIa_Mdl2
Cp->SNIa_Mdu2
Cp->SNIa_bb2
Cp->SNIa_cte2
Cp->SNIa_a2
*/
/***********************************************/
/***********************************************/
/* minimum live time for a given metallicity */
minlivetime
=
star_lifetime
(
StP
[
m
].
Metal
[
NELEMENTS
-
1
],
Cp
->
Mmax
*
SOLAR_MASS
/
All
.
UnitMass_in_g
)
*
All
.
HubbleParam
;
/* maximum live time for a given metallicity */
maxlivetime
=
star_lifetime
(
StP
[
m
].
Metal
[
NELEMENTS
-
1
],
Cp
->
Mmin
*
SOLAR_MASS
/
All
.
UnitMass_in_g
)
*
All
.
HubbleParam
;
//if (P[i].ID==65546)
// printf("(%d) %g the particle 65546 has a max livetime of %g (metal=%g Mmin=%g)\n",ThisTask,All.Time,maxlivetime,StP[m].Metal[NELEMENTS-1],Cp->Mmin);
if
(
All
.
ComovingIntegrationOn
)
{
/* FormationTime on the time line */
Ti0
=
log
(
StP
[
m
].
FormationTime
/
All
.
TimeBegin
)
/
All
.
Timebase_interval
;
/* Beginning of time step on the time line */
Ti1
=
P
[
i
].
Ti_begstep
;
/* End of time step on the time line */
Ti2
=
All
.
Ti_Current
;
#ifdef COSMICTIME
t01
=
get_cosmictime_difference
(
Ti0
,
Ti1
);
t02
=
get_cosmictime_difference
(
Ti0
,
Ti2
);
#endif
}
else
{
t1
=
All
.
TimeBegin
+
(
P
[
i
].
Ti_begstep
*
All
.
Timebase_interval
);
t2
=
All
.
TimeBegin
+
(
All
.
Ti_Current
*
All
.
Timebase_interval
);
t01
=
t1
-
StP
[
m
].
FormationTime
;
t02
=
t2
-
StP
[
m
].
FormationTime
;
}
/* now treat all cases */
do_it
=
1
;
/* beginning of interval */
if
(
t01
>=
minlivetime
)
if
(
t01
>=
maxlivetime
)
do_it
=
0
;
/* nothing to do */
else
m2
=
star_mass_from_age
(
StP
[
m
].
Metal
[
NELEMENTS
-
1
],
t01
/
All
.
HubbleParam
)
*
All
.
HubbleParam
;
else
m2
=
Cp
->
Mmax
*
SOLAR_MASS
/
All
.
UnitMass_in_g
*
All
.
HubbleParam
;
/* end of interval */
if
(
t02
<=
maxlivetime
)
if
(
t02
<=
minlivetime
)
do_it
=
0
;
/* nothing to do */
else
m1
=
star_mass_from_age
(
StP
[
m
].
Metal
[
NELEMENTS
-
1
],
t02
/
All
.
HubbleParam
)
*
All
.
HubbleParam
;
else
m1
=
Cp
->
Mmin
*
SOLAR_MASS
/
All
.
UnitMass_in_g
*
All
.
HubbleParam
;
//printf("Time=%g t01=%g t02=%g id=%d minlivetime=%g maxlivetime=%g \n",All.Time,t01,t02,P[i].ID,minlivetime,maxlivetime);
/* if some of the stars in the SSP explode between t1 and t2 */
if
(
do_it
)
{
Nchimlocal
++
;
StP
[
m
].
Flag
=
1
;
/* mark it as active */
if
(
m1
>
m2
)
{
printf
(
"m1=%g (%g Msol) > m2=%g (%g Msol) !!!
\n\n
"
,
m1
,
m1
*
All
.
UnitMass_in_g
/
SOLAR_MASS
,
m2
,
m2
*
All
.
UnitMass_in_g
/
SOLAR_MASS
);
endrun
(
777002
);
}
M0
=
StP
[
m
].
InitialMass
;
for
(
k
=
0
;
k
<
NELEMENTS
;
k
++
)
metals
[
k
]
=
StP
[
m
].
Metal
[
k
];
/* number of SNIa */
NSNIa
=
SNIa_rate
(
m1
/
All
.
HubbleParam
,
m2
/
All
.
HubbleParam
)
*
M0
/
All
.
HubbleParam
;
/* number of SNII */
NSNII
=
SNII_rate
(
m1
/
All
.
HubbleParam
,
m2
/
All
.
HubbleParam
)
*
M0
/
All
.
HubbleParam
;
NSNIa_totlocal
+=
NSNIa
;
NSNII_totlocal
+=
NSNII
;
/* compute ejectas */
Total_mass_ejection
(
m1
/
All
.
HubbleParam
,
m2
/
All
.
HubbleParam
,
M0
/
All
.
HubbleParam
,
metals
);
StP
[
m
].
TotalEjectedGasMass
=
EjectedMass
[
0
]
*
All
.
HubbleParam
;
/* gas mass */
for
(
k
=
0
;
k
<
NELEMENTS
;
k
++
)
StP
[
m
].
TotalEjectedEltMass
[
k
]
=
EjectedMass
[
k
+
2
]
*
All
.
HubbleParam
;
/* metal mass */
if
(
StP
[
m
].
TotalEjectedGasMass
>
0
)
flag_chimie
=
1
;
/* compute injected energy */
StP
[
m
].
TotalEjectedEgySpec
=
All
.
ChimieSupernovaEnergy
*
(
NSNIa
+
NSNII
)
/
StP
[
m
].
TotalEjectedGasMass
;
EgySNlocal
+=
All
.
ChimieSupernovaEnergy
*
(
NSNIa
+
NSNII
);
/* correct mass particle */
if
(
P
[
i
].
Mass
-
StP
[
m
].
TotalEjectedGasMass
<
0
)
{
printf
(
"mass wants to be less than zero...
\n
"
);
printf
(
"P[i].Mass=%g StP[m].TotalEjectedGasMass=%g
\n
"
,
P
[
i
].
Mass
,
StP
[
m
].
TotalEjectedGasMass
);
endrun
(
777100
);
}
//if (P[i].ID==65546)
// printf("(%d) %g the particle 65546 is here, mass=%g TotalEjectedEltMass=%g m1=%g m2=%g\n",ThisTask,All.Time,P[i].Mass,StP[m].TotalEjectedGasMass,m1,m2);
P
[
i
].
Mass
=
P
[
i
].
Mass
-
StP
[
m
].
TotalEjectedGasMass
;
//float Fe,Mg;
//Fe = StP[m].TotalEjectedEltMass[0];
//Mg = StP[m].TotalEjectedEltMass[1];
}
/******************************************/
/* end do chimie */
/******************************************/
ndone
++
;
if
(
flag_chimie
)
{
for
(
j
=
0
;
j
<
NTask
;
j
++
)
Exportflag
[
j
]
=
0
;
chimie_evaluate
(
i
,
0
);
for
(
j
=
0
;
j
<
NTask
;
j
++
)
{
if
(
Exportflag
[
j
])
{
for
(
k
=
0
;
k
<
3
;
k
++
)
{
ChimieDataIn
[
nexport
].
Pos
[
k
]
=
P
[
i
].
Pos
[
k
];
ChimieDataIn
[
nexport
].
Vel
[
k
]
=
P
[
i
].
Vel
[
k
];
}
ChimieDataIn
[
nexport
].
ID
=
P
[
i
].
ID
;
ChimieDataIn
[
nexport
].
Timestep
=
P
[
i
].
Ti_endstep
-
P
[
i
].
Ti_begstep
;
ChimieDataIn
[
nexport
].
Hsml
=
StP
[
m
].
Hsml
;
ChimieDataIn
[
nexport
].
Density
=
StP
[
m
].
Density
;
ChimieDataIn
[
nexport
].
Volume
=
StP
[
m
].
Volume
;
#ifdef CHIMIE_KINETIC_FEEDBACK
ChimieDataIn
[
nexport
].
NgbMass
=
StP
[
m
].
NgbMass
;
#endif
ChimieDataIn
[
nexport
].
TotalEjectedGasMass
=
StP
[
m
].
TotalEjectedGasMass
;
for
(
k
=
0
;
k
<
NELEMENTS
;
k
++
)
ChimieDataIn
[
nexport
].
TotalEjectedEltMass
[
k
]
=
StP
[
m
].
TotalEjectedEltMass
[
k
];
ChimieDataIn
[
nexport
].
TotalEjectedEgySpec
=
StP
[
m
].
TotalEjectedEgySpec
;
#ifdef WITH_ID_IN_HYDRA
ChimieDataIn
[
nexport
].
ID
=
P
[
i
].
ID
;
#endif
ChimieDataIn
[
nexport
].
Index
=
i
;
ChimieDataIn
[
nexport
].
Task
=
j
;
nexport
++
;
nsend_local
[
j
]
++
;
}
}
}
}
}
}
tend
=
second
();
timecomp
+=
timediff
(
tstart
,
tend
);
qsort
(
ChimieDataIn
,
nexport
,
sizeof
(
struct
chimiedata_in
),
chimie_compare_key
);
for
(
j
=
1
,
noffset
[
0
]
=
0
;
j
<
NTask
;
j
++
)
noffset
[
j
]
=
noffset
[
j
-
1
]
+
nsend_local
[
j
-
1
];
tstart
=
second
();
MPI_Allgather
(
nsend_local
,
NTask
,
MPI_INT
,
nsend
,
NTask
,
MPI_INT
,
MPI_COMM_WORLD
);
tend
=
second
();
timeimbalance
+=
timediff
(
tstart
,
tend
);
/* now do the particles that need to be exported */
for
(
level
=
1
;
level
<
(
1
<<
PTask
);
level
++
)
{
tstart
=
second
();
for
(
j
=
0
;
j
<
NTask
;
j
++
)
nbuffer
[
j
]
=
0
;
for
(
ngrp
=
level
;
ngrp
<
(
1
<<
PTask
);
ngrp
++
)
{
maxfill
=
0
;
for
(
j
=
0
;
j
<
NTask
;
j
++
)
{
if
((
j
^
ngrp
)
<
NTask
)
if
(
maxfill
<
nbuffer
[
j
]
+
nsend
[(
j
^
ngrp
)
*
NTask
+
j
])
maxfill
=
nbuffer
[
j
]
+
nsend
[(
j
^
ngrp
)
*
NTask
+
j
];
}
if
(
maxfill
>=
All
.
BunchSizeChimie
)
break
;
sendTask
=
ThisTask
;
recvTask
=
ThisTask
^
ngrp
;
if
(
recvTask
<
NTask
)
{
if
(
nsend
[
ThisTask
*
NTask
+
recvTask
]
>
0
||
nsend
[
recvTask
*
NTask
+
ThisTask
]
>
0
)
{
/* get the particles */
MPI_Sendrecv
(
&
ChimieDataIn
[
noffset
[
recvTask
]],
nsend_local
[
recvTask
]
*
sizeof
(
struct
chimiedata_in
),
MPI_BYTE
,
recvTask
,
TAG_CHIMIE_A
,
&
ChimieDataGet
[
nbuffer
[
ThisTask
]],
nsend
[
recvTask
*
NTask
+
ThisTask
]
*
sizeof
(
struct
chimiedata_in
),
MPI_BYTE
,
recvTask
,
TAG_CHIMIE_A
,
MPI_COMM_WORLD
,
&
status
);
}
}
for
(
j
=
0
;
j
<
NTask
;
j
++
)
if
((
j
^
ngrp
)
<
NTask
)
nbuffer
[
j
]
+=
nsend
[(
j
^
ngrp
)
*
NTask
+
j
];
}
tend
=
second
();
timecommsumm
+=
timediff
(
tstart
,
tend
);
/* now do the imported particles */
tstart
=
second
();
for
(
j
=
0
;
j
<
nbuffer
[
ThisTask
];
j
++
)
chimie_evaluate
(
j
,
1
);
tend
=
second
();
timecomp
+=
timediff
(
tstart
,
tend
);
/* do a block to measure imbalance */
tstart
=
second
();
MPI_Barrier
(
MPI_COMM_WORLD
);
tend
=
second
();
timeimbalance
+=
timediff
(
tstart
,
tend
);
/* get the result */
tstart
=
second
();
for
(
j
=
0
;
j
<
NTask
;
j
++
)
nbuffer
[
j
]
=
0
;
for
(
ngrp
=
level
;
ngrp
<
(
1
<<
PTask
);
ngrp
++
)
{
maxfill
=
0
;
for
(
j
=
0
;
j
<
NTask
;
j
++
)
{
if
((
j
^
ngrp
)
<
NTask
)
if
(
maxfill
<
nbuffer
[
j
]
+
nsend
[(
j
^
ngrp
)
*
NTask
+
j
])
maxfill
=
nbuffer
[
j
]
+
nsend
[(
j
^
ngrp
)
*
NTask
+
j
];
}
if
(
maxfill
>=
All
.
BunchSizeChimie
)
break
;
sendTask
=
ThisTask
;
recvTask
=
ThisTask
^
ngrp
;
if
(
recvTask
<
NTask
)
{
if
(
nsend
[
ThisTask
*
NTask
+
recvTask
]
>
0
||
nsend
[
recvTask
*
NTask
+
ThisTask
]
>
0
)
{
/* send the results */
MPI_Sendrecv
(
&
ChimieDataResult
[
nbuffer
[
ThisTask
]],
nsend
[
recvTask
*
NTask
+
ThisTask
]
*
sizeof
(
struct
chimiedata_out
),
MPI_BYTE
,
recvTask
,
TAG_CHIMIE_B
,
&
ChimieDataPartialResult
[
noffset
[
recvTask
]],
nsend_local
[
recvTask
]
*
sizeof
(
struct
chimiedata_out
),
MPI_BYTE
,
recvTask
,
TAG_CHIMIE_B
,
MPI_COMM_WORLD
,
&
status
);
/* add the result to the particles */
for
(
j
=
0
;
j
<
nsend_local
[
recvTask
];
j
++
)
{
source
=
j
+
noffset
[
recvTask
];
place
=
ChimieDataIn
[
source
].
Index
;
// for(k = 0; k < 3; k++)
// SphP[place].HydroAccel[k] += HydroDataPartialResult[source].Acc[k];
//
// SphP[place].DtEntropy += HydroDataPartialResult[source].DtEntropy;
//#ifdef FEEDBACK
// SphP[place].DtEgySpecFeedback += HydroDataPartialResult[source].DtEgySpecFeedback;
//#endif
// if(SphP[place].MaxSignalVel < HydroDataPartialResult[source].MaxSignalVel)
// SphP[place].MaxSignalVel = HydroDataPartialResult[source].MaxSignalVel;
//#ifdef COMPUTE_VELOCITY_DISPERSION
// for(k = 0; k < VELOCITY_DISPERSION_SIZE; k++)
// SphP[place].VelocityDispersion[k] += HydroDataPartialResult[source].VelocityDispersion[k];
//#endif
}
}
}
for
(
j
=
0
;
j
<
NTask
;
j
++
)
if
((
j
^
ngrp
)
<
NTask
)
nbuffer
[
j
]
+=
nsend
[(
j
^
ngrp
)
*
NTask
+
j
];
}
tend
=
second
();
timecommsumm
+=
timediff
(
tstart
,
tend
);
level
=
ngrp
-
1
;
}
MPI_Allgather
(
&
ndone
,
1
,
MPI_INT
,
ndonelist
,
1
,
MPI_INT
,
MPI_COMM_WORLD
);
for
(
j
=
0
;
j
<
NTask
;
j
++
)
ntotleft
-=
ndonelist
[
j
];
}
free
(
ndonelist
);
free
(
nsend
);
free
(
nsend_local
);
free
(
nbuffer
);
free
(
noffset
);
/* do final operations on results */
tstart
=
second
();
for
(
i
=
0
;
i
<
N_gas
;
i
++
)
{
if
(
P
[
i
].
Type
==
0
)
{
P
[
i
].
Mass
+=
SphP
[
i
].
dMass
;
SphP
[
i
].
dMass
=
0.
;
}
}
tend
=
second
();
timecomp
+=
timediff
(
tstart
,
tend
);
/* collect some timing information */
MPI_Reduce
(
&
timecomp
,
&
sumt
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Reduce
(
&
timecommsumm
,
&
sumcomm
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Reduce
(
&
timeimbalance
,
&
sumimbalance
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
// if(ThisTask == 0)
// {
// All.CPU_HydCompWalk += sumt / NTask;
// All.CPU_HydCommSumm += sumcomm / NTask;
// All.CPU_HydImbalance += sumimbalance / NTask;
// }
/* collect some chimie informations */
MPI_Reduce
(
&
NSNIa_totlocal
,
&
NSNIa_tot
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Reduce
(
&
NSNII_totlocal
,
&
NSNII_tot
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Reduce
(
&
EgySNlocal
,
&
EgySN
,
1
,
MPI_DOUBLE
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Reduce
(
&
Nchimlocal
,
&
Nchim
,
1
,
MPI_INT
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
#ifdef CHIMIE_THERMAL_FEEDBACK
EgySNThermal
=
EgySN
*
(
1
-
All
.
ChimieKineticFeedbackFraction
);
#else
EgySNThermal
=
0
;
#endif
#ifdef CHIMIE_KINETIC_FEEDBACK
EgySNKinetic
=
EgySN
*
All
.
ChimieKineticFeedbackFraction
;
/* count number of wind particles */
for
(
i
=
0
;
i
<
N_gas
;
i
++
)
{
if
(
P
[
i
].
Type
==
0
)
{
if
(
SphP
[
i
].
WindTime
>=
(
All
.
Time
-
All
.
ChimieWindTime
))
Nwindlocal
++
;
//else
// if (SphP[i].WindTime > All.TimeBegin-2*All.ChimieWindTime)
// Noldwindlocal++;
if
(
SphP
[
i
].
WindFlag
)
Nflaglocal
++
;
}
}
MPI_Reduce
(
&
Nwindlocal
,
&
Nwind
,
1
,
MPI_INT
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Reduce
(
&
Noldwindlocal
,
&
Noldwind
,
1
,
MPI_INT
,
MPI_SUM
,
0
,
MPI_COMM_WORLD
);
MPI_Allreduce
(
&
Nflaglocal
,
&
Nflag
,
1
,
MPI_INT
,
MPI_SUM
,
MPI_COMM_WORLD
);
#else
EgySNKinetic
=
0
;
#endif
/* write some info */
if
(
ThisTask
==
0
)
{
fprintf
(
FdChimie
,
"%15g %10d %15g %15g %15g %15g %15g %10d %10d %10d
\n
"
,
All
.
Time
,
Nchim
,
NSNIa_tot
,
NSNII_tot
,
EgySN
,
EgySNThermal
,
EgySNKinetic
,
Nwind
,
Noldwind
,
Nflag
);
fflush
(
FdChimie
);
}
if
(
Nflag
>
0
)
{
SetMinTimeStepForActives
=
1
;
if
(
ThisTask
==
0
)
fprintf
(
FdLog
,
"%g : !!! set min timestep for active particles !!!
\n
"
,
All
.
Time
);
}
}
/*! This function is the 'core' of the Chemie computation. A target
* particle is specified which may either be local, or reside in the
* communication buffer.
*/
void
chimie_evaluate
(
int
target
,
int
mode
)
{
int
j
,
n
,
startnode
,
numngb
,
numngb_inbox
,
k
;
FLOAT
*
pos
,
*
vel
;
//FLOAT *vel;
//FLOAT mass;
double
h
,
h2
;
double
acc
[
3
];
double
dx
,
dy
,
dz
;
double
wk
,
r
,
r2
,
u
=
0
;
double
hinv
=
1
,
hinv3
;
int
target_stp
;
double
density
;
double
volume
;
#ifdef CHIMIE_KINETIC_FEEDBACK
double
ngbmass
;
double
p
;
#endif
double
aij
;
double
ejectedGasMass
;
double
ejectedEltMass
[
NELEMENTS
];
double
ejectedEgySpec
;
double
mass_k
;
double
NewMass
;
double
fv
,
vi2
,
vj2
;
double
EgySpec
,
NewEgySpec
;
double
DeltaEntropy
;
double
DeltaVel
[
3
];
#ifndef LONGIDS
unsigned
int
id
;
/*!< particle identifier */
#else
unsigned
long
long
id
;
/*!< particle identifier */
#endif
if
(
mode
==
0
)
{
pos
=
P
[
target
].
Pos
;
vel
=
P
[
target
].
Vel
;
id
=
P
[
target
].
ID
;
target_stp
=
P
[
target
].
StPIdx
;
h
=
StP
[
target_stp
].
Hsml
;
density
=
StP
[
target_stp
].
Density
;
volume
=
StP
[
target_stp
].
Volume
;
#ifdef CHIMIE_KINETIC_FEEDBACK
ngbmass
=
StP
[
target_stp
].
NgbMass
;
#endif
ejectedGasMass
=
StP
[
target_stp
].
TotalEjectedGasMass
;
for
(
k
=
0
;
k
<
NELEMENTS
;
k
++
)
ejectedEltMass
[
k
]
=
StP
[
target_stp
].
TotalEjectedEltMass
[
k
];
ejectedEgySpec
=
StP
[
target_stp
].
TotalEjectedEgySpec
;
}
else
{
pos
=
ChimieDataGet
[
target
].
Pos
;
vel
=
ChimieDataGet
[
target
].
Vel
;
id
=
ChimieDataGet
[
target
].
ID
;
h
=
ChimieDataGet
[
target
].
Hsml
;
density
=
ChimieDataGet
[
target
].
Density
;
volume
=
ChimieDataGet
[
target
].
Volume
;
#ifdef CHIMIE_KINETIC_FEEDBACK
ngbmass
=
ChimieDataGet
[
target
].
NgbMass
;
#endif
ejectedGasMass
=
ChimieDataGet
[
target
].
TotalEjectedGasMass
;
for
(
k
=
0
;
k
<
NELEMENTS
;
k
++
)
ejectedEltMass
[
k
]
=
ChimieDataGet
[
target
].
TotalEjectedEltMass
[
k
];
ejectedEgySpec
=
ChimieDataGet
[
target
].
TotalEjectedEgySpec
;
}
/* initialize variables before SPH loop is started */
acc
[
0
]
=
acc
[
1
]
=
acc
[
2
]
=
0
;
vi2
=
0
;
for
(
k
=
0
;
k
<
3
;
k
++
)
vi2
+=
vel
[
k
]
*
vel
[
k
];
h2
=
h
*
h
;
hinv
=
1.0
/
h
;
#ifndef TWODIMS
hinv3
=
hinv
*
hinv
*
hinv
;
#else
hinv3
=
hinv
*
hinv
/
boxSize_Z
;
#endif
/* Now start the actual SPH computation for this particle */
startnode
=
All
.
MaxPart
;
numngb
=
0
;
do
{
numngb_inbox
=
ngb_treefind_variable_for_chimie
(
&
pos
[
0
],
h
,
&
startnode
);
for
(
n
=
0
;
n
<
numngb_inbox
;
n
++
)
{
j
=
Ngblist
[
n
];
dx
=
pos
[
0
]
-
P
[
j
].
Pos
[
0
];
dy
=
pos
[
1
]
-
P
[
j
].
Pos
[
1
];
dz
=
pos
[
2
]
-
P
[
j
].
Pos
[
2
];
#ifdef PERIODIC
/* now find the closest image in the given box size */
if
(
dx
>
boxHalf_X
)
dx
-=
boxSize_X
;
if
(
dx
<
-
boxHalf_X
)
dx
+=
boxSize_X
;
if
(
dy
>
boxHalf_Y
)
dy
-=
boxSize_Y
;
if
(
dy
<
-
boxHalf_Y
)
dy
+=
boxSize_Y
;
if
(
dz
>
boxHalf_Z
)
dz
-=
boxSize_Z
;
if
(
dz
<
-
boxHalf_Z
)
dz
+=
boxSize_Z
;
#endif
r2
=
dx
*
dx
+
dy
*
dy
+
dz
*
dz
;
if
(
r2
<
h2
)
{
numngb
++
;
r
=
sqrt
(
r2
);
u
=
r
*
hinv
;
if
(
u
<
0.5
)
{
wk
=
hinv3
*
(
KERNEL_COEFF_1
+
KERNEL_COEFF_2
*
(
u
-
1
)
*
u
*
u
);
}
else
{
wk
=
hinv3
*
KERNEL_COEFF_5
*
(
1.0
-
u
)
*
(
1.0
-
u
)
*
(
1.0
-
u
);
}
/* normalisation using mass */
aij
=
P
[
j
].
Mass
*
wk
/
density
;
/* normalisation using volume */
/* !!! si on utilise, il faut stoquer une nouvelle variable : OldDensity, car density est modifié plus bas... */
//aij = P[j].Mass/SphP[j].Density*wk/volume;
/* metal injection */
for
(
k
=
0
;
k
<
NELEMENTS
;
k
++
)
{
mass_k
=
SphP
[
j
].
Metal
[
k
]
*
P
[
j
].
Mass
;
/* mass of elt k */
SphP
[
j
].
Metal
[
k
]
=
(
mass_k
+
aij
*
ejectedEltMass
[
k
]
)
/
(
P
[
j
].
Mass
+
aij
*
ejectedGasMass
);
}
/* new mass */
NewMass
=
P
[
j
].
Mass
+
aij
*
ejectedGasMass
;
/* new velocity */
vj2
=
0
;
for
(
k
=
0
;
k
<
3
;
k
++
)
vj2
+=
SphP
[
j
].
VelPred
[
k
]
*
SphP
[
j
].
VelPred
[
k
];
fv
=
sqrt
(
(
P
[
j
].
Mass
/
NewMass
)
+
aij
*
(
ejectedGasMass
/
NewMass
)
*
(
vi2
/
vj2
)
);
for
(
k
=
0
;
k
<
3
;
k
++
)
{
DeltaVel
[
k
]
=
fv
*
SphP
[
j
].
VelPred
[
k
]
-
SphP
[
j
].
VelPred
[
k
];
SphP
[
j
].
VelPred
[
k
]
+=
DeltaVel
[
k
];
P
[
j
].
Vel
[
k
]
+=
DeltaVel
[
k
];
}
/* spec energy at current step */
EgySpec
=
SphP
[
j
].
EntropyPred
/
GAMMA_MINUS1
*
pow
(
SphP
[
j
].
Density
*
a3inv
,
GAMMA_MINUS1
);
/* new egyspec */
NewEgySpec
=
(
EgySpec
)
*
(
P
[
j
].
Mass
/
NewMass
);
/* new density */
SphP
[
j
].
Density
=
SphP
[
j
].
Density
*
NewMass
/
P
[
j
].
Mass
;
/* new entropy */
DeltaEntropy
=
GAMMA_MINUS1
*
NewEgySpec
/
pow
(
SphP
[
j
].
Density
*
a3inv
,
GAMMA_MINUS1
)
-
SphP
[
j
].
EntropyPred
;
SphP
[
j
].
EntropyPred
+=
DeltaEntropy
;
SphP
[
j
].
Entropy
+=
DeltaEntropy
;
#ifdef CHIMIE_THERMAL_FEEDBACK
SphP
[
j
].
DeltaEgySpec
+=
(
1.
-
All
.
ChimieKineticFeedbackFraction
)
*
(
ejectedGasMass
*
ejectedEgySpec
)
*
aij
/
NewMass
;
#endif
#ifdef CHIMIE_KINETIC_FEEDBACK
p
=
(
All
.
ChimieKineticFeedbackFraction
*
ejectedEgySpec
*
ejectedGasMass
)
/
(
0.5
*
ngbmass
*
All
.
ChimieWindSpeed
*
All
.
ChimieWindSpeed
);
double
r
;
r
=
get_Chimie_random_number
(
P
[
j
].
ID
+
id
);
if
(
r
<
p
)
/* we should maybe have a 2d table here... */
{
if
(
SphP
[
j
].
WindTime
<
(
All
.
Time
-
All
.
ChimieWindTime
))
/* not a wind particle */
{
SphP
[
j
].
WindFlag
=
1
;
SphP
[
j
].
WindTime
=
All
.
Time
;
}
}
#endif
#ifdef CHECK_ENTROPY_SIGN
if
((
SphP
[
j
].
EntropyPred
<
0
)
||
(
SphP
[
j
].
Entropy
<
0
))
{
printf
(
"
\n
task=%d: entropy less than zero in chimie_evaluate !
\n
"
,
ThisTask
);
printf
(
"ID=%d Entropy=%g EntropyPred=%g DeltaEntropy=%g
\n
"
,
P
[
j
].
ID
,
SphP
[
j
].
Entropy
,
SphP
[
j
].
EntropyPred
,
DeltaEntropy
);
fflush
(
stdout
);
endrun
(
777003
);
}
#endif
/* store mass diff. */
SphP
[
j
].
dMass
+=
NewMass
-
P
[
j
].
Mass
;
}
}
}
while
(
startnode
>=
0
);
/* Now collect the result at the right place */
if
(
mode
==
0
)
{
// for(k = 0; k < 3; k++)
// SphP[target].HydroAccel[k] = acc[k];
// SphP[target].DtEntropy = dtEntropy;
//#ifdef FEEDBACK
// SphP[target].DtEgySpecFeedback = dtEgySpecFeedback;
//#endif
// SphP[target].MaxSignalVel = maxSignalVel;
//#ifdef COMPUTE_VELOCITY_DISPERSION
// for(k = 0; k < VELOCITY_DISPERSION_SIZE; k++)
// SphP[target].VelocityDispersion[k] = VelocityDispersion[k];
//#endif
}
else
{
// for(k = 0; k < 3; k++)
// HydroDataResult[target].Acc[k] = acc[k];
// HydroDataResult[target].DtEntropy = dtEntropy;
//#ifdef FEEDBACK
// HydroDataResult[target].DtEgySpecFeedback = dtEgySpecFeedback;
//#endif
// HydroDataResult[target].MaxSignalVel = maxSignalVel;
//#ifdef COMPUTE_VELOCITY_DISPERSION
// for(k = 0; k < VELOCITY_DISPERSION_SIZE; k++)
// HydroDataResult[target].VelocityDispersion[k] = VelocityDispersion[k];
//#endif
}
}
/*! This is a comparison kernel for a sort routine, which is used to group
* particles that are going to be exported to the same CPU.
*/
int
chimie_compare_key
(
const
void
*
a
,
const
void
*
b
)
{
if
(((
struct
chimiedata_in
*
)
a
)
->
Task
<
(((
struct
chimiedata_in
*
)
b
)
->
Task
))
return
-
1
;
if
(((
struct
chimiedata_in
*
)
a
)
->
Task
>
(((
struct
chimiedata_in
*
)
b
)
->
Task
))
return
+
1
;
return
0
;
}
#endif
Event Timeline
Log In to Comment