Matching catalogs based on membership (detailed)¶
Here we show the specific steps of matching two catalogs based on proximity
%load_ext autoreload
%autoreload 2
ClCatalogs¶
Given some input data
import numpy as np
from astropy.table import Table
input1 = Table({'ID': ['CL0a', 'CL1a', 'CL2a', 'CL3a', 'CL4a']})
input1['MASS'] = 1e14*np.arange(1, 6)*10
input2 = Table({'ID': ['CL0b', 'CL1b', 'CL2b', 'CL3b']})
input2['MASS'] = 1e14*np.arange(1, 5)*10
display(input1)
display(input2)
input1_mem = Table(
{'ID':[
'MEM0', 'MEM1', 'MEM2', 'MEM3', 'MEM4',
'MEM5', 'MEM6', 'MEM7', 'MEM8', 'MEM9',
'MEM10', 'MEM11', 'MEM12', 'MEM13', 'MEM14'],
'ID_CLUSTER': [
'CL0a', 'CL0a', 'CL0a', 'CL0a', 'CL0a',
'CL1a', 'CL1a', 'CL1a', 'CL1a', 'CL2a',
'CL2a', 'CL2a', 'CL3a', 'CL3a', 'CL4a'],
})
input2_mem = Table(
{'ID':[
'MEM0', 'MEM1', 'MEM2', 'MEM3', 'MEM4',
'MEM5', 'MEM6', 'MEM7', 'MEM8', 'MEM9',
'MEM10', 'MEM11', 'MEM12', 'MEM13'],
'ID_CLUSTER': [
'CL3b', 'CL0b', 'CL0b', 'CL0b', 'CL0b',
'CL1b', 'CL1b', 'CL1b', 'CL1b', 'CL2b',
'CL2b', 'CL2b', 'CL3b', 'CL3b'],
})
input1_mem['RA'] = np.arange(len(input1_mem))*10.0
input2_mem['RA'] = np.arange(len(input2_mem))*10.0
input1_mem['DEC'] = 0.0
input2_mem['DEC'] = 0.0
input1_mem['Z'] = 0.1
input2_mem['Z'] = 0.1
input1_mem['PMEM'] = 1.0
input2_mem['PMEM'] = 1.0
display(input1_mem)
display(input2_mem)
ID | MASS |
---|---|
str4 | float64 |
CL0a | 1000000000000000.0 |
CL1a | 2000000000000000.0 |
CL2a | 3000000000000000.0 |
CL3a | 4000000000000000.0 |
CL4a | 5000000000000000.0 |
ID | MASS |
---|---|
str4 | float64 |
CL0b | 1000000000000000.0 |
CL1b | 2000000000000000.0 |
CL2b | 3000000000000000.0 |
CL3b | 4000000000000000.0 |
ID | ID_CLUSTER | RA | DEC | Z | PMEM |
---|---|---|---|---|---|
str5 | str4 | float64 | float64 | float64 | float64 |
MEM0 | CL0a | 0.0 | 0.0 | 0.1 | 1.0 |
MEM1 | CL0a | 10.0 | 0.0 | 0.1 | 1.0 |
MEM2 | CL0a | 20.0 | 0.0 | 0.1 | 1.0 |
MEM3 | CL0a | 30.0 | 0.0 | 0.1 | 1.0 |
MEM4 | CL0a | 40.0 | 0.0 | 0.1 | 1.0 |
MEM5 | CL1a | 50.0 | 0.0 | 0.1 | 1.0 |
MEM6 | CL1a | 60.0 | 0.0 | 0.1 | 1.0 |
MEM7 | CL1a | 70.0 | 0.0 | 0.1 | 1.0 |
MEM8 | CL1a | 80.0 | 0.0 | 0.1 | 1.0 |
MEM9 | CL2a | 90.0 | 0.0 | 0.1 | 1.0 |
MEM10 | CL2a | 100.0 | 0.0 | 0.1 | 1.0 |
MEM11 | CL2a | 110.0 | 0.0 | 0.1 | 1.0 |
MEM12 | CL3a | 120.0 | 0.0 | 0.1 | 1.0 |
MEM13 | CL3a | 130.0 | 0.0 | 0.1 | 1.0 |
MEM14 | CL4a | 140.0 | 0.0 | 0.1 | 1.0 |
ID | ID_CLUSTER | RA | DEC | Z | PMEM |
---|---|---|---|---|---|
str5 | str4 | float64 | float64 | float64 | float64 |
MEM0 | CL3b | 0.0 | 0.0 | 0.1 | 1.0 |
MEM1 | CL0b | 10.0 | 0.0 | 0.1 | 1.0 |
MEM2 | CL0b | 20.0 | 0.0 | 0.1 | 1.0 |
MEM3 | CL0b | 30.0 | 0.0 | 0.1 | 1.0 |
MEM4 | CL0b | 40.0 | 0.0 | 0.1 | 1.0 |
MEM5 | CL1b | 50.0 | 0.0 | 0.1 | 1.0 |
MEM6 | CL1b | 60.0 | 0.0 | 0.1 | 1.0 |
MEM7 | CL1b | 70.0 | 0.0 | 0.1 | 1.0 |
MEM8 | CL1b | 80.0 | 0.0 | 0.1 | 1.0 |
MEM9 | CL2b | 90.0 | 0.0 | 0.1 | 1.0 |
MEM10 | CL2b | 100.0 | 0.0 | 0.1 | 1.0 |
MEM11 | CL2b | 110.0 | 0.0 | 0.1 | 1.0 |
MEM12 | CL3b | 120.0 | 0.0 | 0.1 | 1.0 |
MEM13 | CL3b | 130.0 | 0.0 | 0.1 | 1.0 |
Create two ClCatalog
objects, they have the same properties of
astropy
tables with additional functionality. For the membership
matching, the main columns to be included are: - id
- must
correspond to id_cluster
in the cluster member catalog. - mass
(or mass proxy) - necessary for proxity matching if
shared_member_fraction
used as preference criteria for unique
matches, default use.
All of the columns can be added when creating the ClCatalog
object
passing them as keys:
cat = ClCatalog('Cat', ra=[0, 1])
or can also be added afterwards:
cat = ClCatalog('Cat')
cat['ra'] = [0, 1]
from clevar.catalog import ClCatalog
c1 = ClCatalog('Cat1', id=input1['ID'], mass=input1['MASS'])
c2 = ClCatalog('Cat2', id=input2['ID'], mass=input2['MASS'])
# Format for nice display
c1['mass'].info.format = '.2e'
c2['mass'].info.format = '.2e'
display(c1)
display(c2)
tags: id(id), mass(mass)
Radius unit: None
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other |
---|---|---|---|---|---|
str4 | float64 | object | object | object | object |
CL0a | 1.00e+15 | None | None | [] | [] |
CL1a | 2.00e+15 | None | None | [] | [] |
CL2a | 3.00e+15 | None | None | [] | [] |
CL3a | 4.00e+15 | None | None | [] | [] |
CL4a | 5.00e+15 | None | None | [] | [] |
tags: id(id), mass(mass)
Radius unit: None
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other |
---|---|---|---|---|---|
str4 | float64 | object | object | object | object |
CL0b | 1.00e+15 | None | None | [] | [] |
CL1b | 2.00e+15 | None | None | [] | [] |
CL2b | 3.00e+15 | None | None | [] | [] |
CL3b | 4.00e+15 | None | None | [] | [] |
The members can be added to the cluster object using the add_members
function. It has a similar instanciating format of a ClCatalog
object, where the columns are added by keyword arguments (the key
id_cluster
is always necessary and must correspond to id
in the
main cluster catalog).
For details see XXX:
from clevar.catalog import MemCatalog
c1.add_members(id=input1_mem['ID'], id_cluster=input1_mem['ID_CLUSTER'],
ra=input1_mem['RA'], dec=input1_mem['DEC'], z=input1_mem['Z'],
pmem=input1_mem['PMEM'])
c2.add_members(id=input2_mem['ID'], id_cluster=input2_mem['ID_CLUSTER'],
ra=input2_mem['RA'], dec=input2_mem['DEC'], z=input2_mem['Z'],
pmem=input2_mem['PMEM'])
display(c1.members)
display(c2.members)
tags: id(id), id_cluster(id_cluster), ra(ra), dec(dec), z(z), pmem(pmem)
id | id_cluster | ra | dec | z | pmem | ind_cl |
---|---|---|---|---|---|---|
str5 | str4 | float64 | float64 | float64 | float64 | int64 |
MEM0 | CL0a | 0.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM1 | CL0a | 10.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM2 | CL0a | 20.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM3 | CL0a | 30.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM4 | CL0a | 40.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM5 | CL1a | 50.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM6 | CL1a | 60.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM7 | CL1a | 70.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM8 | CL1a | 80.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM9 | CL2a | 90.0 | 0.0 | 0.1 | 1.0 | 2 |
MEM10 | CL2a | 100.0 | 0.0 | 0.1 | 1.0 | 2 |
MEM11 | CL2a | 110.0 | 0.0 | 0.1 | 1.0 | 2 |
MEM12 | CL3a | 120.0 | 0.0 | 0.1 | 1.0 | 3 |
MEM13 | CL3a | 130.0 | 0.0 | 0.1 | 1.0 | 3 |
MEM14 | CL4a | 140.0 | 0.0 | 0.1 | 1.0 | 4 |
tags: id(id), id_cluster(id_cluster), ra(ra), dec(dec), z(z), pmem(pmem)
id | id_cluster | ra | dec | z | pmem | ind_cl |
---|---|---|---|---|---|---|
str5 | str4 | float64 | float64 | float64 | float64 | int64 |
MEM0 | CL3b | 0.0 | 0.0 | 0.1 | 1.0 | 3 |
MEM1 | CL0b | 10.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM2 | CL0b | 20.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM3 | CL0b | 30.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM4 | CL0b | 40.0 | 0.0 | 0.1 | 1.0 | 0 |
MEM5 | CL1b | 50.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM6 | CL1b | 60.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM7 | CL1b | 70.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM8 | CL1b | 80.0 | 0.0 | 0.1 | 1.0 | 1 |
MEM9 | CL2b | 90.0 | 0.0 | 0.1 | 1.0 | 2 |
MEM10 | CL2b | 100.0 | 0.0 | 0.1 | 1.0 | 2 |
MEM11 | CL2b | 110.0 | 0.0 | 0.1 | 1.0 | 2 |
MEM12 | CL3b | 120.0 | 0.0 | 0.1 | 1.0 | 3 |
MEM13 | CL3b | 130.0 | 0.0 | 0.1 | 1.0 | 3 |
The catalogs can also be read directly from files, for more details see catalogs.ipynb.
Matching¶
Import the MembershipMatch
and create a object for matching
from clevar.match import MembershipMatch
mt = MembershipMatch()
Prepare the matching object¶
Before matching the clusters it is necessary to match the members catalogs and then filling the clusters with information about of the shared members.
The matching of members can be done by id
if both member catalogs
share the same id
s or by angular proximity.
The first step is to prepare each catalog with the matching configuration:
delta_z
: Defines redshift window for matching. The possible values are:'cat'
: uses redshift properties of the catalog'spline.filename'
: interpolates data in'filename'
assuming (z, zmin, zmax) formatfloat
: usesdelta_z*(1+z)
None
: does not use z
match_radius
: Radius of the catalog to be used in the matching. If'cat'
uses the radius in the catalog, else must be in format'value unit'
. (ex:'1 arcsec'
,'1 Mpc'
)
In this case, because one of the configuraion radius has physical units,
we also need a cosmology (cosmo
) object to convert it to angular
size (this is done internally).
To match the members by id
, just run the function:
%%time
mt.match_members(c1.members, c2.members, method='id')
28 members were matched.
CPU times: user 404 µs, sys: 634 µs, total: 1.04 ms
Wall time: 1.03 ms
To match the members by angular proximity you also have to provide: -
radius
(str
, None
). Radius for matching, with format
'value unit'
(ex: 1 arcsec
, 1 Mpc
). -
cosmo
(clevar.Cosmology
, None
). Cosmology object for when
radius has physical units. Then call the same function with these
arguments
from clevar.cosmology import AstroPyCosmology
mt.match_members(c1.members, c2.members, method='angular_distance',
radius='0.1 kpc', cosmo=AstroPyCosmology())
## ClCatalog 1
## Prep mt_cols
* zmin|zmax set to -1|10
* ang radius from set scale
## ClCatalog 2
## Prep mt_cols
* zmin|zmax set to -1|10
* ang radius from set scale
## Multiple match (catalog 1)
Finding candidates (members)
* 14/15 objects matched.
## Multiple match (catalog 2)
Finding candidates (members)
* 14/14 objects matched.
## Finding unique matches of catalog 1
Unique Matches (members)
* 14/15 objects matched.
## Finding unique matches of catalog 2
Unique Matches (members)
* 14/14 objects matched.
Cross Matches (members)
* 14/15 objects matched.
Cross Matches (members)
* 14/14 objects matched.
28 members were matched.
This function adds a matched_mems
attribute to your matching object
(mt.matched_mems
in this case) that contains the indices of the
matched members. This attribute can be saved and loaded so you don’t
have to redo this step. Just use the functions:
mt.save_matched_members(filename='mem_mt.txt', overwrite=False)
mt.load_matched_members(filename='mem_mt.txt')
Now we fill the catalogs with the information regarding the matched
members. In this step, each cluster catalog will have a ClData
table
in its mt_input
attibute with the number of members in each cluster
(nmem
) and a dictionary containing the number of shaded objects with
the clusters of the other catalog (shared_mems
).
If pmem
is provided to the members, these quantities are computed
as:
\(nmem=\sum_{cluster\;members} Pmem_i\)
\(shared\_mems=\sum_{shared\;members} Pmem_i\)
mt.fill_shared_members(c1, c2)
display(c1.mt_input)
display(c2.mt_input)
share_mems | nmem |
---|---|
object | float64 |
{'CL3b': 1.0, 'CL0b': 4.0} | 5.0 |
{'CL1b': 4.0} | 4.0 |
{'CL2b': 3.0} | 3.0 |
{'CL3b': 2.0} | 2.0 |
{} | 1.0 |
share_mems | nmem |
---|---|
object | float64 |
{'CL0a': 4.0} | 4.0 |
{'CL1a': 4.0} | 4.0 |
{'CL2a': 3.0} | 3.0 |
{'CL3a': 2.0, 'CL0a': 1.0} | 3.0 |
Again, these results can be saved and loaded so you don’t have to redo this step. Just use the functions:
mt.save_shared_members(c1, c2, fileprefix='mem_share')
mt.load_shared_members(c1, c2, fileprefix='mem_share')
Once this step is done, you can actually start matching the clusters.
Multiple matching¶
The next step is to match the catalogs and store all candidates that pass the matching criteria.
%%time
mt.multiple(c1, c2)
mt.multiple(c2, c1)
Finding candidates (Cat1)
* 4/5 objects matched.
Finding candidates (Cat2)
* 4/4 objects matched.
CPU times: user 1.98 ms, sys: 0 ns, total: 1.98 ms
Wall time: 1.86 ms
This will fill the mt_multi_self
and mt_multi_other
columns:
display(c1)
display(c2)
tags: id(id), mass(mass)
Radius unit: None
mt_input | |||||||
---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | share_mems | nmem |
str4 | float64 | object | object | object | object | object | float64 |
CL0a | 1.00e+15 | None | None | ['CL0b', 'CL3b'] | ['CL0b', 'CL3b'] | {'CL3b': 1.0, 'CL0b': 4.0} | 5.0 |
CL1a | 2.00e+15 | None | None | ['CL1b'] | ['CL1b'] | {'CL1b': 4.0} | 4.0 |
CL2a | 3.00e+15 | None | None | ['CL2b'] | ['CL2b'] | {'CL2b': 3.0} | 3.0 |
CL3a | 4.00e+15 | None | None | ['CL3b'] | ['CL3b'] | {'CL3b': 2.0} | 2.0 |
CL4a | 5.00e+15 | None | None | [] | [] | {} | 1.0 |
tags: id(id), mass(mass)
Radius unit: None
mt_input | |||||||
---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | share_mems | nmem |
str4 | float64 | object | object | object | object | object | float64 |
CL0b | 1.00e+15 | None | None | ['CL0a'] | ['CL0a'] | {'CL0a': 4.0} | 4.0 |
CL1b | 2.00e+15 | None | None | ['CL1a'] | ['CL1a'] | {'CL1a': 4.0} | 4.0 |
CL2b | 3.00e+15 | None | None | ['CL2a'] | ['CL2a'] | {'CL2a': 3.0} | 3.0 |
CL3b | 4.00e+15 | None | None | ['CL3a', 'CL0a'] | ['CL3a', 'CL0a'] | {'CL3a': 2.0, 'CL0a': 1.0} | 3.0 |
Unique matching¶
Once all candidates are stored in each catalog, we can find the best
candidates. You can also pass the argument: - preference
: In cases
where there are multiple matched, how the best candidate will be chosen.
Options are: 'more_massive'
, 'angular_proximity'
,
'redshift_proximity'
, 'shared_member_fraction'
(default value).
%%time
mt.unique(c1, c2, preference='shared_member_fraction')
mt.unique(c2, c1, preference='shared_member_fraction')
Unique Matches (Cat1)
* 4/5 objects matched.
Unique Matches (Cat2)
* 4/4 objects matched.
CPU times: user 695 µs, sys: 919 µs, total: 1.61 ms
Wall time: 1.56 ms
This will fill the matching columns: - mt_self
: Best candidate found
- mt_other
: Best candidate found by the other catalog -
mt_frac_self
: Fraction of shared members with the best candidate
found - mt_frac_other
: Fraction of shared members by the best
candidate found by the other catalog, relative to the other catalog
If pmem
is present in the members catalogs, the shared fractions are
computed by:
\(\frac{\sum_{shared\;members}Pmem_i}{\sum_{cluster\;members}Pmem_i}\)
\(\frac{\sum_{shared\;members}Pmem_i}{\sum_{cluster\;members}Pmem_i}\)
display(c1)
display(c2)
tags: id(id), mass(mass)
Radius unit: None
mt_input | |||||||||
---|---|---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_frac_self | mt_frac_other | share_mems | nmem |
str4 | float64 | object | object | object | object | float64 | float64 | object | float64 |
CL0a | 1.00e+15 | CL0b | CL0b | ['CL0b', 'CL3b'] | ['CL0b', 'CL3b'] | 0.8 | 1.0 | {'CL3b': 1.0, 'CL0b': 4.0} | 5.0 |
CL1a | 2.00e+15 | CL1b | CL1b | ['CL1b'] | ['CL1b'] | 1.0 | 1.0 | {'CL1b': 4.0} | 4.0 |
CL2a | 3.00e+15 | CL2b | CL2b | ['CL2b'] | ['CL2b'] | 1.0 | 1.0 | {'CL2b': 3.0} | 3.0 |
CL3a | 4.00e+15 | CL3b | CL3b | ['CL3b'] | ['CL3b'] | 1.0 | 0.6666666666666666 | {'CL3b': 2.0} | 2.0 |
CL4a | 5.00e+15 | None | None | [] | [] | 0.0 | 0.0 | {} | 1.0 |
tags: id(id), mass(mass)
Radius unit: None
mt_input | |||||||||
---|---|---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_frac_other | mt_frac_self | share_mems | nmem |
str4 | float64 | object | object | object | object | float64 | float64 | object | float64 |
CL0b | 1.00e+15 | CL0a | CL0a | ['CL0a'] | ['CL0a'] | 0.8 | 1.0 | {'CL0a': 4.0} | 4.0 |
CL1b | 2.00e+15 | CL1a | CL1a | ['CL1a'] | ['CL1a'] | 1.0 | 1.0 | {'CL1a': 4.0} | 4.0 |
CL2b | 3.00e+15 | CL2a | CL2a | ['CL2a'] | ['CL2a'] | 1.0 | 1.0 | {'CL2a': 3.0} | 3.0 |
CL3b | 4.00e+15 | CL3a | CL3a | ['CL3a', 'CL0a'] | ['CL3a', 'CL0a'] | 1.0 | 0.6666666666666666 | {'CL3a': 2.0, 'CL0a': 1.0} | 3.0 |
Cross matching¶
If you want to make sure the same pair was found in both directions:
c1.cross_match()
c2.cross_match()
This will fill the mt_cross
column:
display(c1)
display(c2)
tags: id(id), mass(mass)
Radius unit: None
mt_input | ||||||||||
---|---|---|---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_frac_self | mt_frac_other | mt_cross | share_mems | nmem |
str4 | float64 | object | object | object | object | float64 | float64 | object | object | float64 |
CL0a | 1.00e+15 | CL0b | CL0b | ['CL0b', 'CL3b'] | ['CL0b', 'CL3b'] | 0.8 | 1.0 | CL0b | {'CL3b': 1.0, 'CL0b': 4.0} | 5.0 |
CL1a | 2.00e+15 | CL1b | CL1b | ['CL1b'] | ['CL1b'] | 1.0 | 1.0 | CL1b | {'CL1b': 4.0} | 4.0 |
CL2a | 3.00e+15 | CL2b | CL2b | ['CL2b'] | ['CL2b'] | 1.0 | 1.0 | CL2b | {'CL2b': 3.0} | 3.0 |
CL3a | 4.00e+15 | CL3b | CL3b | ['CL3b'] | ['CL3b'] | 1.0 | 0.6666666666666666 | CL3b | {'CL3b': 2.0} | 2.0 |
CL4a | 5.00e+15 | None | None | [] | [] | 0.0 | 0.0 | None | {} | 1.0 |
tags: id(id), mass(mass)
Radius unit: None
mt_input | ||||||||||
---|---|---|---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_frac_other | mt_frac_self | mt_cross | share_mems | nmem |
str4 | float64 | object | object | object | object | float64 | float64 | object | object | float64 |
CL0b | 1.00e+15 | CL0a | CL0a | ['CL0a'] | ['CL0a'] | 0.8 | 1.0 | CL0a | {'CL0a': 4.0} | 4.0 |
CL1b | 2.00e+15 | CL1a | CL1a | ['CL1a'] | ['CL1a'] | 1.0 | 1.0 | CL1a | {'CL1a': 4.0} | 4.0 |
CL2b | 3.00e+15 | CL2a | CL2a | ['CL2a'] | ['CL2a'] | 1.0 | 1.0 | CL2a | {'CL2a': 3.0} | 3.0 |
CL3b | 4.00e+15 | CL3a | CL3a | ['CL3a', 'CL0a'] | ['CL3a', 'CL0a'] | 1.0 | 0.6666666666666666 | CL3a | {'CL3a': 2.0, 'CL0a': 1.0} | 3.0 |
Save and Load¶
The results of the matching can easily be saved and load using
ClEvaR
tools:
mt.save_matches(c1, c2, out_dir='temp', overwrite=True)
mt.load_matches(c1, c2, out_dir='temp')
display(c1)
display(c2)
Cat1
<< ClEvar used in matching: 0.13.2 >>
* Total objects: 5
* multiple (self): 4
* multiple (other): 4
* unique (self): 4
* unique (other): 4
* cross: 4
Cat2
<< ClEvar used in matching: 0.13.2 >>
* Total objects: 4
* multiple (self): 4
* multiple (other): 4
* unique (self): 4
* unique (other): 4
* cross: 4
tags: id(id), mass(mass)
Radius unit: None
mt_input | ||||||||||
---|---|---|---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_frac_self | mt_frac_other | mt_cross | share_mems | nmem |
str4 | float64 | object | object | object | object | float64 | float64 | object | object | float64 |
CL0a | 1.00e+15 | CL0b | CL0b | ['CL0b', 'CL3b'] | ['CL0b', 'CL3b'] | 0.8 | 1.0 | CL0b | {'CL3b': 1.0, 'CL0b': 4.0} | 5.0 |
CL1a | 2.00e+15 | CL1b | CL1b | ['CL1b'] | ['CL1b'] | 1.0 | 1.0 | CL1b | {'CL1b': 4.0} | 4.0 |
CL2a | 3.00e+15 | CL2b | CL2b | ['CL2b'] | ['CL2b'] | 1.0 | 1.0 | CL2b | {'CL2b': 3.0} | 3.0 |
CL3a | 4.00e+15 | CL3b | CL3b | ['CL3b'] | ['CL3b'] | 1.0 | 0.6666666666666666 | CL3b | {'CL3b': 2.0} | 2.0 |
CL4a | 5.00e+15 | None | None | [] | [] | 0.0 | 0.0 | None | {} | 1.0 |
tags: id(id), mass(mass)
Radius unit: None
mt_input | ||||||||||
---|---|---|---|---|---|---|---|---|---|---|
id | mass | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_frac_other | mt_frac_self | mt_cross | share_mems | nmem |
str4 | float64 | object | object | object | object | float64 | float64 | object | object | float64 |
CL0b | 1.00e+15 | CL0a | CL0a | ['CL0a'] | ['CL0a'] | 0.8 | 1.0 | CL0a | {'CL0a': 4.0} | 4.0 |
CL1b | 2.00e+15 | CL1a | CL1a | ['CL1a'] | ['CL1a'] | 1.0 | 1.0 | CL1a | {'CL1a': 4.0} | 4.0 |
CL2b | 3.00e+15 | CL2a | CL2a | ['CL2a'] | ['CL2a'] | 1.0 | 1.0 | CL2a | {'CL2a': 3.0} | 3.0 |
CL3b | 4.00e+15 | CL3a | CL3a | ['CL3a', 'CL0a'] | ['CL3a', 'CL0a'] | 1.0 | 0.6666666666666666 | CL3a | {'CL3a': 2.0, 'CL0a': 1.0} | 3.0 |
Getting Matched Pairs¶
There is functionality inbuilt in clevar
to plot some results of the
matching, such as: - Recovery rates - Distances (anguar and redshift) of
cluster centers - Scaling relations (mass, redshift, …) for those cases,
check the match_metrics.ipynb and match_metrics_advanced.ipynb
notebooks.
If those do not provide your needs, you can get directly the matched pairs of clusters:
from clevar.match import get_matched_pairs
mt1, mt2 = get_matched_pairs(c1, c2, 'cross')
These will be catalogs with the corresponding matched pairs:
import pylab as plt
plt.scatter(mt1['mass'], mt2['mass'])
<matplotlib.collections.PathCollection at 0x7f02dbe02730>

Members of matched pairs¶
The members also carry the information on the matched clusters. The
column match
shows to which clusters of the other catalog this
member also belongs. The column in_mt_sample
says if those clusters
are presented in the matched sample:
mt1.members
tags: id(id), id_cluster(id_cluster), ra(ra), dec(dec), z(z), pmem(pmem)
id | id_cluster | ra | dec | z | pmem | ind_cl | match | mt_self | mt_other | mt_multi_self | mt_multi_other | mt_cross | in_mt_sample |
---|---|---|---|---|---|---|---|---|---|---|---|---|---|
str5 | str4 | float64 | float64 | float64 | float64 | int64 | object | object | object | object | object | object | bool |
MEM0 | CL0a | 0.0 | 0.0 | 0.1 | 1.0 | 0 | ['CL3b'] | MEM0 | MEM0 | ['MEM0'] | ['MEM0'] | MEM0 | True |
MEM1 | CL0a | 10.0 | 0.0 | 0.1 | 1.0 | 0 | ['CL0b'] | MEM1 | MEM1 | ['MEM1'] | ['MEM1'] | MEM1 | True |
MEM2 | CL0a | 20.0 | 0.0 | 0.1 | 1.0 | 0 | ['CL0b'] | MEM2 | MEM2 | ['MEM2'] | ['MEM2'] | MEM2 | True |
MEM3 | CL0a | 30.0 | 0.0 | 0.1 | 1.0 | 0 | ['CL0b'] | MEM3 | MEM3 | ['MEM3'] | ['MEM3'] | MEM3 | True |
MEM4 | CL0a | 40.0 | 0.0 | 0.1 | 1.0 | 0 | ['CL0b'] | MEM4 | MEM4 | ['MEM4'] | ['MEM4'] | MEM4 | True |
MEM5 | CL1a | 50.0 | 0.0 | 0.1 | 1.0 | 1 | ['CL1b'] | MEM5 | MEM5 | ['MEM5'] | ['MEM5'] | MEM5 | True |
MEM6 | CL1a | 60.0 | 0.0 | 0.1 | 1.0 | 1 | ['CL1b'] | MEM6 | MEM6 | ['MEM6'] | ['MEM6'] | MEM6 | True |
MEM7 | CL1a | 70.0 | 0.0 | 0.1 | 1.0 | 1 | ['CL1b'] | MEM7 | MEM7 | ['MEM7'] | ['MEM7'] | MEM7 | True |
MEM8 | CL1a | 80.0 | 0.0 | 0.1 | 1.0 | 1 | ['CL1b'] | MEM8 | MEM8 | ['MEM8'] | ['MEM8'] | MEM8 | True |
MEM9 | CL2a | 90.0 | 0.0 | 0.1 | 1.0 | 2 | ['CL2b'] | MEM9 | MEM9 | ['MEM9'] | ['MEM9'] | MEM9 | True |
MEM10 | CL2a | 100.0 | 0.0 | 0.1 | 1.0 | 2 | ['CL2b'] | MEM10 | MEM10 | ['MEM10'] | ['MEM10'] | MEM10 | True |
MEM11 | CL2a | 110.0 | 0.0 | 0.1 | 1.0 | 2 | ['CL2b'] | MEM11 | MEM11 | ['MEM11'] | ['MEM11'] | MEM11 | True |
MEM12 | CL3a | 120.0 | 0.0 | 0.1 | 1.0 | 3 | ['CL3b'] | MEM12 | MEM12 | ['MEM12'] | ['MEM12'] | MEM12 | True |
MEM13 | CL3a | 130.0 | 0.0 | 0.1 | 1.0 | 3 | ['CL3b'] | MEM13 | MEM13 | ['MEM13'] | ['MEM13'] | MEM13 | True |
Outputing matched catalogs¶
To save the current catalogs, you can use the write
inbuilt
function:
c1.write('c1_temp.fits', overwrite=True)
This will allow you to save the catalog with its current labels and matching information.
Outputing matching information to original catalogs¶
Assuming your input data came from initial files, clevar
also
provides functions create output files that combine all the information
on them with the matching results.
To add the matching information to an input catalog, use:
from clevar.match import output_catalog_with_matching
output_catalog_with_matching('input_catalog.fits', 'output_catalog.fits', c1)
note:
input_catalog.fits
must have the same number of rows thatc1
.
To create a matched catalog containig all columns of both input catalogs, use:
from clevar.match import output_matched_catalog
output_matched_catalog('input_catalog1.fits', 'input_catalog2.fits',
'output_catalog.fits', c1, c2, matching_type='cross')
where matching_type
must be cross
, cat1
or cat2
.
note:
input_catalog1.fits
must have the same number of rows thatc1
(and the same forc2
).