@@ -111,7 +111,7 @@ template <class MASKType, class B0Type, class DKTType, class PredictType>
111
111
class Processor { MEMALIGN(Processor)
112
112
public:
113
113
Processor (const Eigen::MatrixXd& b, const bool ols, const int iter,
114
- MASKType* mask_image, B0Type* b0_image, DKTType* dkt_image, PredictType* predict_image) :
114
+ const MASKType& mask_image, const B0Type& b0_image, const DKTType& dkt_image, const PredictType& predict_image) :
115
115
mask_image (mask_image),
116
116
b0_image (b0_image),
117
117
dkt_image (dkt_image),
@@ -128,9 +128,9 @@ class Processor { MEMALIGN(Processor)
128
128
template <class DWIType , class DTType >
129
129
void operator () (DWIType& dwi_image, DTType& dt_image)
130
130
{
131
- if (mask_image) {
132
- assign_pos_of (dwi_image, 0 , 3 ).to (* mask_image);
133
- if (!mask_image-> value ())
131
+ if (mask_image. valid () ) {
132
+ assign_pos_of (dwi_image, 0 , 3 ).to (mask_image);
133
+ if (!mask_image. value ())
134
134
return ;
135
135
}
136
136
@@ -157,34 +157,34 @@ class Processor { MEMALIGN(Processor)
157
157
dt_image.value () = p[dt_image.index (3 )];
158
158
}
159
159
160
- if (b0_image) {
161
- assign_pos_of (dwi_image, 0 , 3 ).to (* b0_image);
162
- b0_image-> value () = exp (p[6 ]);
160
+ if (b0_image. valid () ) {
161
+ assign_pos_of (dwi_image, 0 , 3 ).to (b0_image);
162
+ b0_image. value () = exp (p[6 ]);
163
163
}
164
164
165
- if (dkt_image) {
166
- assign_pos_of (dwi_image, 0 , 3 ).to (* dkt_image);
165
+ if (dkt_image. valid () ) {
166
+ assign_pos_of (dwi_image, 0 , 3 ).to (dkt_image);
167
167
double adc_sq = (p[0 ]+p[1 ]+p[2 ])*(p[0 ]+p[1 ]+p[2 ])/9.0 ;
168
- for (auto l = Loop (3 )(* dkt_image); l; ++l) {
169
- dkt_image-> value () = p[dkt_image-> index (3 )+7 ]/adc_sq;
168
+ for (auto l = Loop (3 )(dkt_image); l; ++l) {
169
+ dkt_image. value () = p[dkt_image. index (3 )+7 ]/adc_sq;
170
170
}
171
171
}
172
172
173
- if (predict_image) {
174
- assign_pos_of (dwi_image, 0 , 3 ).to (* predict_image);
173
+ if (predict_image. valid () ) {
174
+ assign_pos_of (dwi_image, 0 , 3 ).to (predict_image);
175
175
dwi = (b*p).array ().exp ();
176
- for (auto l = Loop (3 )(* predict_image); l; ++l) {
177
- predict_image-> value () = dwi[predict_image-> index (3 )];
176
+ for (auto l = Loop (3 )(predict_image); l; ++l) {
177
+ predict_image. value () = dwi[predict_image. index (3 )];
178
178
}
179
179
}
180
180
181
181
}
182
182
183
183
private:
184
- copy_ptr< MASKType> mask_image;
185
- copy_ptr< B0Type> b0_image;
186
- copy_ptr< DKTType> dkt_image;
187
- copy_ptr< PredictType> predict_image;
184
+ MASKType mask_image;
185
+ B0Type b0_image;
186
+ DKTType dkt_image;
187
+ PredictType predict_image;
188
188
Eigen::VectorXd dwi;
189
189
Eigen::VectorXd p;
190
190
Eigen::VectorXd w;
@@ -196,7 +196,7 @@ class Processor { MEMALIGN(Processor)
196
196
};
197
197
198
198
template <class MASKType , class B0Type , class DKTType , class PredictType >
199
- inline Processor<MASKType, B0Type, DKTType, PredictType> processor (const Eigen::MatrixXd& b, const bool ols, const int iter, MASKType* mask_image, B0Type* b0_image, DKTType* dkt_image, PredictType* predict_image) {
199
+ inline Processor<MASKType, B0Type, DKTType, PredictType> processor (const Eigen::MatrixXd& b, const bool ols, const int iter, const MASKType& mask_image, const B0Type& b0_image, const DKTType& dkt_image, const PredictType& predict_image) {
200
200
return { b, ols, iter, mask_image, b0_image, dkt_image, predict_image };
201
201
}
202
202
@@ -205,11 +205,11 @@ void run ()
205
205
auto dwi = Header::open (argument[0 ]).get_image <value_type>();
206
206
auto grad = DWI::get_DW_scheme (dwi);
207
207
208
- Image<bool >* mask = nullptr ;
208
+ Image<bool > mask;
209
209
auto opt = get_options (" mask" );
210
210
if (opt.size ()) {
211
- mask = new Image<bool > (Image< bool > ::open (opt[0 ][0 ]) );
212
- check_dimensions (dwi, * mask, 0 , 3 );
211
+ mask = Image<bool >::open (opt[0 ][0 ]);
212
+ check_dimensions (dwi, mask, 0 , 3 );
213
213
}
214
214
215
215
bool ols = get_options (" ols" ).size ();
@@ -223,27 +223,27 @@ void run ()
223
223
DWI::stash_DW_scheme (header, grad);
224
224
PhaseEncoding::clear_scheme (header);
225
225
226
- Image<value_type>* predict = nullptr ;
226
+ Image<value_type> predict;
227
227
opt = get_options (" predicted_signal" );
228
228
if (opt.size ())
229
- predict = new Image<value_type> (Image<value_type> ::create (opt[0 ][0 ], header) );
229
+ predict = Image<value_type>::create (opt[0 ][0 ], header);
230
230
231
231
header.size (3 ) = 6 ;
232
232
auto dt = Image<value_type>::create (argument[1 ], header);
233
233
234
- Image<value_type>* b0 = nullptr ;
234
+ Image<value_type> b0;
235
235
opt = get_options (" b0" );
236
236
if (opt.size ()) {
237
237
header.ndim () = 3 ;
238
- b0 = new Image<value_type> (Image<value_type> ::create (opt[0 ][0 ], header) );
238
+ b0 = Image<value_type>::create (opt[0 ][0 ], header);
239
239
}
240
240
241
- Image<value_type>* dkt = nullptr ;
241
+ Image<value_type> dkt;
242
242
opt = get_options (" dkt" );
243
243
if (opt.size ()) {
244
244
header.ndim () = 4 ;
245
245
header.size (3 ) = 15 ;
246
- dkt = new Image<value_type> (Image<value_type> ::create (opt[0 ][0 ], header) );
246
+ dkt = Image<value_type>::create (opt[0 ][0 ], header);
247
247
}
248
248
249
249
Eigen::MatrixXd b = -DWI::grad2bmatrix<double > (grad, opt.size ()>0 );
0 commit comments