mirror of
https://github.com/gwm17/catima.git
synced 2024-11-23 02:38:51 -05:00
commit
f212a7409b
35
catima.cpp
35
catima.cpp
|
@ -87,6 +87,23 @@ double dedx_from_range(Projectile &p, double T, const Material &t, const Config
|
|||
return p.A/range_spline.derivative(T);
|
||||
}
|
||||
|
||||
std::vector<double> dedx_from_range(Projectile &p, const std::vector<double> &T, const Material &t, const Config &c){
|
||||
auto data = _storage.Get(p,t,c);
|
||||
Interpolator range_spline(energy_table.values,data.range.data(),energy_table.num);
|
||||
|
||||
std::vector<double> dedx;
|
||||
dedx.reserve(T.size());
|
||||
for(auto e:T){
|
||||
if(e<catima::Ezero){
|
||||
dedx.push_back(0.0);
|
||||
}
|
||||
else{
|
||||
dedx.push_back(p.A/range_spline.derivative(e));
|
||||
}
|
||||
}
|
||||
return dedx;
|
||||
}
|
||||
|
||||
double range_straggling(Projectile &p, double T, const Material &t, const Config &c){
|
||||
double r=0;
|
||||
auto data = _storage.Get(p,t,c);
|
||||
|
@ -161,6 +178,24 @@ double energy_out(Projectile &p, double T, const Material &t, const Config &c){
|
|||
return energy_out(T,t.thickness(),range_spline);
|
||||
}
|
||||
|
||||
std::vector<double> energy_out(Projectile &p, const std::vector<double> &T, const Material &t, const Config &c){
|
||||
auto data = _storage.Get(p,t,c);
|
||||
Interpolator range_spline(energy_table.values,data.range.data(),energy_table.num);
|
||||
|
||||
std::vector<double> eout;
|
||||
eout.reserve(T.size());
|
||||
for(auto e:T){
|
||||
if(e<catima::Ezero){
|
||||
eout.push_back(0.0);
|
||||
}
|
||||
else{
|
||||
eout.push_back(energy_out(e,t.thickness(),range_spline));
|
||||
}
|
||||
}
|
||||
|
||||
return eout;
|
||||
}
|
||||
|
||||
Result calculate(Projectile &p, const Material &t, const Config &c){
|
||||
Result res;
|
||||
double T = p.T;
|
||||
|
|
19
catima.h
19
catima.h
|
@ -18,6 +18,7 @@
|
|||
#define CPPATIMA_H
|
||||
|
||||
#include <utility>
|
||||
#include <vector>
|
||||
|
||||
// #define NDEBUG
|
||||
#include "catima/build_config.h"
|
||||
|
@ -70,6 +71,15 @@ namespace catima{
|
|||
*/
|
||||
double dedx_from_range(Projectile &p, double T, const Material &t, const Config &c=default_config);
|
||||
|
||||
/**
|
||||
* returns the dEdx calculated from range spline as derivative
|
||||
* @param p - Projectile
|
||||
* @param T - energy vector
|
||||
* @param mat - Material
|
||||
* @return range
|
||||
*/
|
||||
std::vector<double> dedx_from_range(Projectile &p, const std::vector<double> &T, const Material &t, const Config &c=default_config);
|
||||
|
||||
/**
|
||||
* returns the range straggling of the Projectile in Material from spline
|
||||
* @param p - Projectile
|
||||
|
@ -144,6 +154,15 @@ namespace catima{
|
|||
*/
|
||||
double energy_out(Projectile &p, double T, const Material &t, const Config &c=default_config);
|
||||
|
||||
/**
|
||||
* calculates outcoming energy
|
||||
* @p - Projectile
|
||||
* @t - Material
|
||||
* @param T - incoming energy vector
|
||||
* @return outcoming energy after the material in Mev/u
|
||||
*/
|
||||
std::vector<double> energy_out(Projectile &p, const std::vector<double> &T, const Material &t, const Config &c=default_config);
|
||||
|
||||
/**
|
||||
* calculates all observables for projectile passing material
|
||||
* @param p - Projectile
|
||||
|
|
21
catima.pyx
21
catima.pyx
|
@ -7,6 +7,7 @@
|
|||
"""
|
||||
|
||||
cimport catimac
|
||||
from libcpp.vector cimport vector
|
||||
from enum import IntEnum
|
||||
import numpy
|
||||
|
||||
|
@ -391,13 +392,13 @@ def projectile_range(Projectile projectile, Material material, energy = None, Co
|
|||
|
||||
def dedx_from_range(Projectile projectile, Material material, energy = None, Config config = default_config):
|
||||
if(isinstance(energy,numpy.ndarray)):
|
||||
res = numpy.empty(energy.size)
|
||||
for i,e in enumerate(energy):
|
||||
res[i] = catimac.dedx_from_range(projectile.cbase, e, material.cbase, config.cbase)
|
||||
return res
|
||||
res = catimac.dedx_from_range(projectile.cbase, <vector[double]>energy.tolist(), material.cbase, config.cbase)
|
||||
return numpy.asarray(res);
|
||||
if(isinstance(energy,list)):
|
||||
return catimac.dedx_from_range(projectile.cbase, <vector[double]>energy, material.cbase, config.cbase)
|
||||
if(energy is None):
|
||||
energy = projectile.T()
|
||||
return catimac.dedx_from_range(projectile.cbase, energy, material.cbase, config.cbase);
|
||||
return catimac.dedx_from_range(projectile.cbase, <double>energy, material.cbase, config.cbase);
|
||||
|
||||
def domega2de(Projectile projectile, Material material, energy = None, Config config = default_config):
|
||||
if(isinstance(energy,numpy.ndarray)):
|
||||
|
@ -441,13 +442,13 @@ def domega2dx(Projectile projectile, Material material, energy = None, Config co
|
|||
|
||||
def energy_out(Projectile projectile, Material material, energy = None, Config config = default_config):
|
||||
if(isinstance(energy,numpy.ndarray)):
|
||||
res = numpy.empty(energy.size)
|
||||
for i,e in enumerate(energy):
|
||||
res[i] = catimac.energy_out(projectile.cbase, e, material.cbase, config.cbase)
|
||||
return res
|
||||
res = catimac.energy_out(projectile.cbase, <vector[double]>energy.tolist(), material.cbase, config.cbase)
|
||||
return numpy.asarray(res)
|
||||
if(isinstance(energy,list)):
|
||||
return catimac.energy_out(projectile.cbase, <vector[double]>energy, material.cbase, config.cbase)
|
||||
if(energy is None):
|
||||
energy = projectile.T()
|
||||
return catimac.energy_out(projectile.cbase, energy, material.cbase, config.cbase)
|
||||
return catimac.energy_out(projectile.cbase, <double>energy, material.cbase, config.cbase)
|
||||
|
||||
def sezi_dedx_e(Projectile projectile, Target t):
|
||||
return catimac.sezi_dedx_e(projectile.cbase, t.cbase)
|
||||
|
|
|
@ -73,7 +73,9 @@ cdef extern from "catima/catima.h" namespace "catima":
|
|||
cdef double domega2dx(Projectile &p, double T, const Material &mat, const Config &c)
|
||||
cdef double range(Projectile &p, double T, const Material &t, const Config &c);
|
||||
cdef double dedx_from_range(Projectile &p, double T, const Material &t, const Config &c);
|
||||
cdef vector[double] dedx_from_range(Projectile &p, vector[double] &T, const Material &t, const Config &c);
|
||||
cdef double energy_out(Projectile &p, double T, const Material &t, const Config &c);
|
||||
cdef vector[double] energy_out(Projectile &p, vector[double] &T, const Material &t, const Config &c);
|
||||
|
||||
cdef double domega2de(Projectile &p, double T, const Material &t, const Config &c);
|
||||
cdef double da2de(Projectile &p, double T, const Material &t, const Config &c);
|
||||
|
|
|
@ -144,8 +144,6 @@ class TestStructures(unittest.TestCase):
|
|||
|
||||
def test_material_calculation(self):
|
||||
water = catima.get_material(catima.material.WATER)
|
||||
water2 = catima.get_material(catima.material.WATER)
|
||||
water2.I(78.0)
|
||||
p = catima.Projectile(1,1)
|
||||
|
||||
p(1000)
|
||||
|
@ -154,15 +152,10 @@ class TestStructures(unittest.TestCase):
|
|||
self.assertAlmostEqual(res.dEdxi,2.23,1)
|
||||
self.assertAlmostEqual(res["dEdxi"],2.23,1)
|
||||
self.assertAlmostEqual(res.dEdxi,res2,3)
|
||||
res3 = catima.calculate(p,water2)
|
||||
self.assertTrue(res.dEdxi>res3.dEdxi)
|
||||
|
||||
res = catima.calculate(p(500),water)
|
||||
res2 = catima.dedx(p,water)
|
||||
self.assertAlmostEqual(res.dEdxi,2.76,1)
|
||||
self.assertAlmostEqual(res.dEdxi,res2,3)
|
||||
res3 = catima.calculate(p,water2)
|
||||
self.assertTrue(res.dEdxi>res3.dEdxi)
|
||||
|
||||
res = catima.calculate(p(9),water)
|
||||
res2 = catima.dedx(p,water)
|
||||
|
@ -187,7 +180,31 @@ class TestStructures(unittest.TestCase):
|
|||
self.assertAlmostEqual(res["Eout"],997.077,1)
|
||||
self.assertAlmostEqual(res.Eout,res2,3)
|
||||
|
||||
def test_layer_calculateion(self):
|
||||
def test_eout_list(self):
|
||||
graphite = catima.get_material(6)
|
||||
graphite.thickness(0.5)
|
||||
p = catima.Projectile(12,6)
|
||||
energies = [100,500,1000]
|
||||
res = catima.calculate(p(1000),graphite)
|
||||
self.assertAlmostEqual(res.Eout,997.077,1)
|
||||
res2 = catima.energy_out(p,graphite,energy=energies)
|
||||
self.assertEqual(len(res2),len(energies))
|
||||
self.assertAlmostEqual(res2[2], 997.077,1)
|
||||
self.assertAlmostEqual(res2[0], catima.calculate(p(energies[0]),graphite).Eout ,1)
|
||||
self.assertAlmostEqual(res2[1], catima.calculate(p(energies[1]),graphite).Eout ,1)
|
||||
|
||||
def test_dedx_from_range_list(self):
|
||||
graphite = catima.get_material(6)
|
||||
graphite.thickness(0.5)
|
||||
p = catima.Projectile(12,6)
|
||||
energies = [100,500,1000]
|
||||
res2 = catima.dedx_from_range(p,graphite,energy=energies)
|
||||
self.assertEqual(len(res2),len(energies))
|
||||
for i,e in enumerate(energies):
|
||||
r = catima.dedx_from_range(p, graphite, energy=e)
|
||||
self.assertAlmostEqual(res2[i], r, 0.1)
|
||||
|
||||
def test_layer_calculation(self):
|
||||
p = catima.Projectile(12,6)
|
||||
water = catima.get_material(catima.material.WATER)
|
||||
water.thickness(10.0)
|
||||
|
@ -261,9 +278,8 @@ class TestStructures(unittest.TestCase):
|
|||
graphite = catima.get_material(6)
|
||||
graphite.thickness(1.0)
|
||||
data = catima.get_data(p, water)
|
||||
self.assertEqual(catima.max_storage_data,50) # assuming 50, this has to be changed manually
|
||||
self.assertEqual(catima.max_storage_data,100) # assuming 50, this has to be changed manually
|
||||
r = catima.storage_info()
|
||||
print(r)
|
||||
|
||||
#self.assertAlmostEqual(catima.da2de(p,water,et[100]),data[2][100],6)
|
||||
#self.assertAlmostEqual(catima.da2de(p,water,et[400]),data[2][400],6)
|
||||
|
|
|
@ -358,6 +358,34 @@ const lest::test specification[] =
|
|||
c.z_effective = z_eff_type::atima14;
|
||||
EXPECT(z_eff_atima14(92,1900,13) == z_effective(p_u(1900.),t,c));
|
||||
#endif
|
||||
},
|
||||
CASE("vector_inputs"){
|
||||
catima::Projectile p{12,6,6,1000};
|
||||
catima::Material water({
|
||||
{1.00794,1,2},
|
||||
{15.9994,8,1}
|
||||
});
|
||||
catima::Material graphite;
|
||||
graphite.add_element(12,6,1);
|
||||
graphite.density(2.0);
|
||||
graphite.thickness(0.5);
|
||||
|
||||
auto res = catima::calculate(p,graphite);
|
||||
EXPECT( res.Eout == approx(997.07,01));
|
||||
|
||||
std::vector<double> energies{100,500,1000};
|
||||
auto res2 = catima::energy_out(p,energies, graphite);
|
||||
EXPECT(res2.size()==energies.size());
|
||||
EXPECT(res2[2] == approx(997.07,01));
|
||||
EXPECT(res2[0] == approx(catima::energy_out(p,energies[0],graphite),0.1));
|
||||
EXPECT(res2[1] == approx(catima::energy_out(p,energies[1],graphite),0.1));
|
||||
EXPECT(res2[2] == approx(catima::energy_out(p,energies[2],graphite),0.1));
|
||||
|
||||
auto res3 = catima::dedx_from_range(p,energies,graphite);
|
||||
EXPECT(res3.size()==energies.size());
|
||||
EXPECT(res3[0] == approx(catima::dedx_from_range(p,energies[0],graphite),0.1));
|
||||
EXPECT(res3[1] == approx(catima::dedx_from_range(p,energies[1],graphite),0.1));
|
||||
EXPECT(res3[2] == approx(catima::dedx_from_range(p,energies[2],graphite),0.1));
|
||||
}
|
||||
|
||||
};
|
||||
|
|
Loading…
Reference in New Issue
Block a user