1 2 3 4 5 6 7 8 9 10 11 12 13 14 15 16 17 18 19 20 21 22 23 24 25 26 27 28 29 30 31 32 33 34 35 36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58 59 60 61 62 63 64 65 66 67 68 69 70 71 72 73 74 75 76 77 78 79 80 81 82 83 84 85 86 87 88 89 90 91 92 93 94 95 96 97 98 99 100 101 102 103 104 105 106 107 108 109 110 111 112 113 114 115 116 117 118 119 120 121 122 123 124 125 126 127 128 129 130 131 132 133 134 135 136 137 138 139 140 141 142 143 144 145 146 147 148 149 150 151 152
|
/* Copyright (c) 2008-2025 the MRtrix3 contributors.
*
* This Source Code Form is subject to the terms of the Mozilla Public
* License, v. 2.0. If a copy of the MPL was not distributed with this
* file, You can obtain one at http://mozilla.org/MPL/2.0/.
*
* Covered Software is provided under this License on an "as is"
* basis, without warranty of any kind, either expressed, implied, or
* statutory, including, without limitation, warranties that the
* Covered Software is free of defects, merchantable, fit for a
* particular purpose or non-infringing.
* See the Mozilla Public License v. 2.0 for more details.
*
* For more details, see http://www.mrtrix.org/.
*/
#include "command.h"
#include "exception.h"
#include "image.h"
#include "mrtrix.h"
#include "progressbar.h"
#include "algo/loop.h"
#include "dwi/gradient.h"
#include "dwi/shells.h"
#include "math/math.h"
#include "math/SH.h"
#include "math/ZSH.h"
using namespace MR;
using namespace App;
void usage ()
{
AUTHOR = "J-Donald Tournier (jdtournier@gmail.com)";
SYNOPSIS = "Generate an appropriate response function from the image data for spherical deconvolution";
DESCRIPTION
+ Math::SH::encoding_description;
ARGUMENTS
+ Argument ("SH", "the spherical harmonic decomposition of the diffusion-weighted images").type_image_in()
+ Argument ("mask", "the mask containing the voxels from which to estimate the response function").type_image_in()
+ Argument ("directions", "a 4D image containing the direction vectors along which to estimate the response function").type_image_in()
+ Argument ("response", "the output axially-symmetric spherical harmonic coefficients").type_file_out();
OPTIONS
+ Option ("lmax", "specify the maximum harmonic degree of the response function to estimate")
+ Argument ("value").type_integer (0, 20)
+ Option ("dump", "dump the m=0 SH coefficients from all voxels in the mask to the output file, rather than their mean")
+ Argument ("file").type_file_out();
}
using value_type = double;
void run ()
{
auto SH = Image<value_type>::open(argument[0]);
Math::SH::check (SH);
auto mask = Image<bool>::open(argument[1]);
auto dir = Image<value_type>::open(argument[2]).with_direct_io();
int lmax = get_option_value ("lmax", Math::SH::LforN (SH.size(3)));
check_dimensions (SH, mask, 0, 3);
check_dimensions (SH, dir, 0, 3);
if (dir.ndim() != 4)
throw Exception ("input direction image \"" + std::string (argument[2]) + "\" must be a 4D image");
if (dir.size(3) != 3)
throw Exception ("input direction image \"" + std::string (argument[2]) + "\" must contain precisely 3 volumes");
Eigen::VectorXd delta;
Eigen::VectorXd response = Eigen::VectorXd::Zero (Math::ZSH::NforL (lmax));
size_t count = 0;
File::OFStream dump_stream;
auto opt = get_options ("dump");
if (opt.size())
dump_stream.open (opt[0][0]);
Eigen::Matrix<value_type,Eigen::Dynamic,1,0,64> AL (lmax+1);
Math::Legendre::Plm_sph (AL, lmax, 0, value_type (1.0));
auto loop = Loop ("estimating response function", SH, 0, 3);
for (auto l = loop(mask, SH, dir); l; ++l) {
if (!mask.value())
continue;
Eigen::Vector3d d = dir.row(3);
if (!d.allFinite()) {
WARN ("voxel with invalid direction [ " + str(dir.index(0)) + " " + str(dir.index(1)) + " " + str(dir.index(2)) + " ]; skipping");
continue;
}
d.normalize();
// Uncertainty regarding Eigen's behaviour when normalizing a zero vector; may change behaviour between versions
if (!d.allFinite() || !d.squaredNorm()) {
WARN ("voxel with zero direction [ " + str(dir.index(0)) + " " + str(dir.index(1)) + " " + str(dir.index(2)) + " ]; skipping");
continue;
}
Math::SH::delta (delta, d, lmax);
for (int l = 0; l <= lmax; l += 2) {
value_type d_dot_s = 0.0;
value_type d_dot_d = 0.0;
for (int m = -l; m <= l; ++m) {
size_t i = Math::SH::index (l,m);
SH.index(3) = i;
value_type s = SH.value();
// TODO: currently this does NOT handle the non-orthonormal basis
d_dot_s += s*delta[i];
d_dot_d += Math::pow2 (delta[i]);
}
value_type val = AL[l] * d_dot_s / d_dot_d;
response[Math::ZSH::index(l)] += val;
if (dump_stream.is_open())
dump_stream << val << " ";
}
if (dump_stream.is_open())
dump_stream << "\n";
++count;
}
response /= count;
if (std::string(argument[3]) == "-") {
for (ssize_t n = 0; n < response.size(); ++n)
std::cout << response[n] << " ";
std::cout << "\n";
}
else {
save_vector(response, argument[3]);
}
}
|