浏览代码

Working with all at once kernel calculation. Much faster.

master
T-bone 7 年前
父节点
当前提交
0c5bbd8c51
共有 3 个文件被更改,包括 58 次插入85 次删除
  1. 4
    4
      examples/KernelV0.cpp
  2. 7
    10
      include/KernelV0.h
  3. 47
    71
      src/KernelV0.cpp

+ 4
- 4
examples/KernelV0.cpp 查看文件

34
 
34
 
35
     // Transmitter loops
35
     // Transmitter loops
36
     auto Tx1 = CircularLoop(21, 15, 100, 100);
36
     auto Tx1 = CircularLoop(21, 15, 100, 100);
37
-    auto Tx2 = CircularLoop(21, 15, 100, 120);
37
+    auto Tx2 = CircularLoop(21, 15, 100, 124.8);
38
     //auto Tx1 = CircularLoop(60, 15, 0, 0); // was 60
38
     //auto Tx1 = CircularLoop(60, 15, 0, 0); // was 60
39
 
39
 
40
     auto Kern = KernelV0::NewSP();
40
     auto Kern = KernelV0::NewSP();
45
 
45
 
46
         Kern->SetIntegrationSize( (Vector3r() << 200,200,200).finished() );
46
         Kern->SetIntegrationSize( (Vector3r() << 200,200,200).finished() );
47
         Kern->SetIntegrationOrigin( (Vector3r() << 0,0,0).finished() );
47
         Kern->SetIntegrationOrigin( (Vector3r() << 0,0,0).finished() );
48
-        Kern->SetTolerance( 1e-12 );
48
+        Kern->SetTolerance( 1e-10 );
49
         //Kern->SetTolerance( .55 ) ; // 1%
49
         //Kern->SetTolerance( .55 ) ; // 1%
50
 
50
 
51
         Kern->SetPulseDuration(0.020);
51
         Kern->SetPulseDuration(0.020);
63
 
63
 
64
     // We could, I suppose, take the earth model in here? For non-linear that
64
     // We could, I suppose, take the earth model in here? For non-linear that
65
     // may be more natural to work with?
65
     // may be more natural to work with?
66
-    std::vector<std::string> tx = {std::string("Coil 1")};
66
+    std::vector<std::string> tx = {std::string("Coil 1"), std::string("Coil 2") };
67
     std::vector<std::string> rx = {std::string("Coil 1")};
67
     std::vector<std::string> rx = {std::string("Coil 1")};
68
-    Kern->CalculateK0( tx, rx, false );
68
+    Kern->CalculateK0( tx, rx, true );
69
 
69
 
70
 }
70
 }
71
 
71
 

+ 7
- 10
include/KernelV0.h 查看文件

10
 /**
10
 /**
11
  * @file
11
  * @file
12
  * @date      11/11/2016 01:47:34 PM
12
  * @date      11/11/2016 01:47:34 PM
13
- * @version   $Id$
14
  * @author    Trevor Irons (ti)
13
  * @author    Trevor Irons (ti)
15
  * @email     tirons@egi.utah.edu
14
  * @email     tirons@egi.utah.edu
16
  * @copyright Copyright (c) 2016, University of Utah
15
  * @copyright Copyright (c) 2016, University of Utah
17
  * @copyright Copyright (c) 2016, Lemma Software, LLC
16
  * @copyright Copyright (c) 2016, Lemma Software, LLC
17
+ * @copyright Copyright (c) 2008, Colorado School of Mines
18
  */
18
  */
19
 #ifndef  KERNELV0_INC
19
 #ifndef  KERNELV0_INC
20
 #define  KERNELV0_INC
20
 #define  KERNELV0_INC
225
          */
225
          */
226
         VectorXcr f( const Vector3r& r, const Real& volume , const Vector3cr& Ht, const Vector3cr& Hr);
226
         VectorXcr f( const Vector3r& r, const Real& volume , const Vector3cr& Ht, const Vector3cr& Hr);
227
 
227
 
228
-//         Complex ComputeV0Cell(const Vector3cr& Bt, const Vector3cr& Br, const Real& vol,
229
-//             const Real& phi);
230
-
231
-        Complex ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
232
-                const Real& sintheta, const Real& phase, const Real& Mn0Abs,
233
-                const Real& vol);
228
+//         Complex ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
229
+//                 const Real& sintheta, const Real& phase, const Real& Mn0Abs,
230
+//                 const Real& vol);
234
 
231
 
235
         EllipticB EllipticFieldRep (const Vector3cr& B, const Vector3r& B0hat);
232
         EllipticB EllipticFieldRep (const Vector3cr& B, const Vector3r& B0hat);
236
 
233
 
237
         Vector3r ComputeMn0(const Real& Porosity, const Vector3r& B0);
234
         Vector3r ComputeMn0(const Real& Porosity, const Vector3r& B0);
238
 
235
 
239
-        Complex IntegrateOnOctreeGrid( const int& iq, bool vtkOutput=false );
236
+        void IntegrateOnOctreeGrid( bool vtkOutput=false );
240
 
237
 
241
         /**
238
         /**
242
          *  Recursive call to integrate a function on an adaptive Octree Grid.
239
          *  Recursive call to integrate a function on an adaptive Octree Grid.
257
          *  HyperOctree, which is useful for visualization.
254
          *  HyperOctree, which is useful for visualization.
258
          */
255
          */
259
         void EvaluateKids2(  const Vector3r& size, const int& level, const Vector3r& cpos,
256
         void EvaluateKids2(  const Vector3r& size, const int& level, const Vector3r& cpos,
260
-                            const Complex& parentVal, vtkHyperOctree* octree, vtkHyperOctreeCursor* curse );
257
+                            const VectorXcr& parentVal, vtkHyperOctree* octree, vtkHyperOctreeCursor* curse );
261
 
258
 
262
         void GetPosition( vtkHyperOctreeCursor* Cursor, Real* p );
259
         void GetPosition( vtkHyperOctreeCursor* Cursor, Real* p );
263
         #endif
260
         #endif
293
         std::map< std::string , std::shared_ptr< EMEarth1D > >             EMEarths;
290
         std::map< std::string , std::shared_ptr< EMEarth1D > >             EMEarths;
294
 
291
 
295
         #ifdef LEMMAUSEVTK
292
         #ifdef LEMMAUSEVTK
296
-        std::map< int, Complex  >                 LeafDict;
293
+        std::map< int, VectorXcr  >               LeafDict;
297
         std::map< int, int     >                  LeafDictIdx;
294
         std::map< int, int     >                  LeafDictIdx;
298
         std::map< int, Real     >                 LeafDictErr;
295
         std::map< int, Real     >                 LeafDictErr;
299
         #endif
296
         #endif

+ 47
- 71
src/KernelV0.cpp 查看文件

10
 /**
10
 /**
11
  * @file
11
  * @file
12
  * @date      11/11/2016 01:47:25 PM
12
  * @date      11/11/2016 01:47:25 PM
13
- * @version   $Id$
14
  * @author    Trevor Irons (ti)
13
  * @author    Trevor Irons (ti)
15
  * @email     tirons@egi.utah.edu
14
  * @email     tirons@egi.utah.edu
16
  * @copyright Copyright (c) 2016, University of Utah
15
  * @copyright Copyright (c) 2016, University of Utah
17
  * @copyright Copyright (c) 2016, Lemma Software, LLC
16
  * @copyright Copyright (c) 2016, Lemma Software, LLC
17
+ * @copyright Copyright (c) 2008, Colorado School of Mines
18
  */
18
  */
19
 
19
 
20
 
20
 
145
             std::cout << "Layer " << ilay << std::endl; //<< " q " << iq << std::endl;
145
             std::cout << "Layer " << ilay << std::endl; //<< " q " << iq << std::endl;
146
             Size(2) = Interfaces(ilay+1) - Interfaces(ilay);
146
             Size(2) = Interfaces(ilay+1) - Interfaces(ilay);
147
             Origin(2) = Interfaces(ilay);
147
             Origin(2) = Interfaces(ilay);
148
-            //for (int iq=0; iq< PulseI.size(); ++iq) {
149
-            //    Ip = PulseI(iq);
150
-                //Kern(ilay, iq) =
151
-            IntegrateOnOctreeGrid( 0, vtkOutput );
152
-            //}
148
+            IntegrateOnOctreeGrid( vtkOutput );
153
         }
149
         }
154
         std::cout << "\rFinished KERNEL\n";
150
         std::cout << "\rFinished KERNEL\n";
155
-        std::cout << "real\n";
151
+        std::cout << "#real\n";
156
         std::cout << Kern.real() << std::endl;
152
         std::cout << Kern.real() << std::endl;
157
-        std::cout << "imag\n";
153
+        std::cout << "#imag\n";
158
         std::cout << Kern.imag() << std::endl;
154
         std::cout << Kern.imag() << std::endl;
159
     }
155
     }
160
 
156
 
162
     //       Class:  KernelV0
158
     //       Class:  KernelV0
163
     //      Method:  IntegrateOnOctreeGrid
159
     //      Method:  IntegrateOnOctreeGrid
164
     //--------------------------------------------------------------------------------------
160
     //--------------------------------------------------------------------------------------
165
-    Complex KernelV0::IntegrateOnOctreeGrid( const int& iq, bool vtkOutput) {
161
+    void KernelV0::IntegrateOnOctreeGrid( bool vtkOutput) {
166
 
162
 
167
         Vector3r cpos = Origin + Size/2.;
163
         Vector3r cpos = Origin + Size/2.;
168
 
164
 
179
                 oct->SetSize( Size(0), Size(1), Size(2) );
175
                 oct->SetSize( Size(0), Size(1), Size(2) );
180
             vtkHyperOctreeCursor* curse = oct->NewCellCursor();
176
             vtkHyperOctreeCursor* curse = oct->NewCellCursor();
181
                 curse->ToRoot();
177
                 curse->ToRoot();
182
-            EvaluateKids2( Size, 0, cpos, 1e6, oct, curse );
178
+            EvaluateKids2( Size, 0, cpos, VectorXcr::Ones(PulseI.size()), oct, curse );
179
+
180
+            for (int iq=0; iq<PulseI.size(); ++iq) {
183
 
181
 
184
             // Fill in leaf data
182
             // Fill in leaf data
185
             vtkDoubleArray* kr = vtkDoubleArray::New();
183
             vtkDoubleArray* kr = vtkDoubleArray::New();
204
 
202
 
205
             //Real LeafVol(0);
203
             //Real LeafVol(0);
206
             for (auto leaf : LeafDict) {
204
             for (auto leaf : LeafDict) {
207
-                kr->InsertTuple1( leaf.first, std::real(leaf.second) );
208
-                ki->InsertTuple1( leaf.first, std::imag(leaf.second) );
209
-                km->InsertTuple1( leaf.first, std::abs(leaf.second) );
205
+                kr->InsertTuple1( leaf.first, std::real(leaf.second(iq)) );
206
+                ki->InsertTuple1( leaf.first, std::imag(leaf.second(iq)) );
207
+                km->InsertTuple1( leaf.first, std::abs(leaf.second(iq)) );
210
                 kid->InsertTuple1( leaf.first, leaf.first );
208
                 kid->InsertTuple1( leaf.first, leaf.first );
211
                 //LeafVol += std::real(leaf.second);
209
                 //LeafVol += std::real(leaf.second);
212
             }
210
             }
216
                 kerr->InsertTuple1( leaf.first, leaf.second );
214
                 kerr->InsertTuple1( leaf.first, leaf.second );
217
             }
215
             }
218
 
216
 
219
-            oct->GetLeafData()->AddArray(kr);
220
-            oct->GetLeafData()->AddArray(ki);
221
-            oct->GetLeafData()->AddArray(km);
222
-            oct->GetLeafData()->AddArray(kid);
223
-            oct->GetLeafData()->AddArray(kerr);
217
+            auto kri = oct->GetLeafData()->AddArray(kr);
218
+            auto kii = oct->GetLeafData()->AddArray(ki);
219
+            auto kmi = oct->GetLeafData()->AddArray(km);
220
+            auto kidi = oct->GetLeafData()->AddArray(kid);
221
+            auto keri = oct->GetLeafData()->AddArray(kerr);
224
 
222
 
225
             auto write = vtkXMLHyperOctreeWriter::New();
223
             auto write = vtkXMLHyperOctreeWriter::New();
226
                 //write.SetDataModeToAscii()
224
                 //write.SetDataModeToAscii()
231
                 write->Write();
229
                 write->Write();
232
                 write->Delete();
230
                 write->Delete();
233
 
231
 
234
-            //kerr->Delete();
232
+            oct->GetLeafData()->RemoveArray( kri );
233
+            oct->GetLeafData()->RemoveArray( kii );
234
+            oct->GetLeafData()->RemoveArray( kmi );
235
+            oct->GetLeafData()->RemoveArray( kidi );
236
+            oct->GetLeafData()->RemoveArray( keri );
237
+
238
+            kerr->Delete();
235
             kid->Delete();
239
             kid->Delete();
236
             kr->Delete();
240
             kr->Delete();
237
             ki->Delete();
241
             ki->Delete();
238
             km->Delete();
242
             km->Delete();
243
+
244
+            }
245
+
239
             curse->Delete();
246
             curse->Delete();
240
             oct->Delete();
247
             oct->Delete();
241
         #else
248
         #else
245
         }
252
         }
246
         std::cout << "\nVOLSUM=" << VOLSUM << "\tActual=" <<  Size(0)*Size(1)*Size(2)
253
         std::cout << "\nVOLSUM=" << VOLSUM << "\tActual=" <<  Size(0)*Size(1)*Size(2)
247
                   << "\tDifference=" << VOLSUM - (Size(0)*Size(1)*Size(2)) <<  std::endl;
254
                   << "\tDifference=" << VOLSUM - (Size(0)*Size(1)*Size(2)) <<  std::endl;
248
-        std::cout << "nleaves\t" << nleaves << std::endl;
249
-        std::cout << "KSUM\t" << SUM << std::endl;
250
-
251
-        return SUM;
252
     }
255
     }
253
 
256
 
254
     //--------------------------------------------------------------------------------------
257
     //--------------------------------------------------------------------------------------
270
         Real Mn0Abs = Mn0.norm();
273
         Real Mn0Abs = Mn0.norm();
271
 
274
 
272
         // Compute phase delay
275
         // Compute phase delay
273
-        // TODO add transmiiter phase and delay time phase!
274
-        Real phase = EBR.zeta+EBT.zeta;
276
+        // TODO add transmiiter current phase and delay induced apparent time phase!
277
+        Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) + (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
278
+        Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
275
 
279
 
276
         // Calcuate vector of all responses
280
         // Calcuate vector of all responses
277
         VectorXcr F = VectorXcr::Zero( PulseI.size() );
281
         VectorXcr F = VectorXcr::Zero( PulseI.size() );
278
         for (int iq=0; iq<PulseI.size(); ++iq) {
282
         for (int iq=0; iq<PulseI.size(); ++iq) {
279
             // Compute the tipping angle
283
             // Compute the tipping angle
280
             Real sintheta = std::sin(0.5*GAMMA*PulseI(iq)*Taup*std::abs(EBT.alpha-EBT.beta));
284
             Real sintheta = std::sin(0.5*GAMMA*PulseI(iq)*Taup*std::abs(EBT.alpha-EBT.beta));
281
-            F(iq) = ComputeV0Cell(EBT, EBR, sintheta, phase, Mn0Abs, volume);
285
+            F(iq) = -volume*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
282
         }
286
         }
283
         return F;
287
         return F;
284
     }
288
     }
287
 //     //       Class:  KernelV0
291
 //     //       Class:  KernelV0
288
 //     //      Method:  ComputeV0Cell
292
 //     //      Method:  ComputeV0Cell
289
 //     //--------------------------------------------------------------------------------------
293
 //     //--------------------------------------------------------------------------------------
290
-//     Complex KernelV0::ComputeV0Cell(const Vector3cr& Bt,
291
-//                 const Vector3cr& Br, const Real& vol, const Real& phi) {
292
-//
293
-//         // Compute the elliptic fields
294
+//     Complex KernelV0::ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
295
+//                 const Real& sintheta, const Real& phase, const Real& Mn0Abs,
296
+//                 const Real& vol) {
297
+//         // earth response of receiver adjoint field
294
 //         Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
298
 //         Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
295
-//         Vector3r B0 = SigmaModel->GetMagneticField();
296
-//
297
-//         // Elliptic representation
298
-//         EllipticB EBT = EllipticFieldRep(Bt, B0hat);
299
-//         EllipticB EBR = EllipticFieldRep(Br, B0hat);
300
-//
301
-//         // Compute Mn0
302
-//         Vector3r Mn0 = ComputeMn0(phi, B0);
303
-//         Real Mn0Abs = Mn0.norm();
304
-//
305
-//         // Compute phase delay, TODO add transmiiter phase and delay time phase!
306
-//         Real phase = EBR.zeta+EBT.zeta;
307
-//
308
-//         Real sintheta = std::sin(0.5*GAMMA*Ip*Taup*std::abs(EBT.alpha-EBT.beta));
309
-//         return 0; ComputeV0Cell(EBT, EBR, sintheta, phase, Mn0Abs, vol);
299
+//         Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
300
+//         Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) + (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
301
+//         return -vol*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
310
 //     }
302
 //     }
311
 
303
 
312
     //--------------------------------------------------------------------------------------
304
     //--------------------------------------------------------------------------------------
313
     //       Class:  KernelV0
305
     //       Class:  KernelV0
314
     //      Method:  ComputeV0Cell
306
     //      Method:  ComputeV0Cell
315
     //--------------------------------------------------------------------------------------
307
     //--------------------------------------------------------------------------------------
316
-    Complex KernelV0::ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
317
-                const Real& sintheta, const Real& phase, const Real& Mn0Abs,
318
-                const Real& vol) {
319
-        // earth response of receiver adjoint field
320
-        Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
321
-        Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
322
-        Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) + (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
323
-        return -vol*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
324
-    }
325
-
326
-    //--------------------------------------------------------------------------------------
327
-    //       Class:  KernelV0
328
-    //      Method:  ComputeV0Cell
329
-    //--------------------------------------------------------------------------------------
330
     Vector3r KernelV0::ComputeMn0(const Real& Porosity, const Vector3r& B0) {
308
     Vector3r KernelV0::ComputeMn0(const Real& Porosity, const Vector3r& B0) {
331
         Real chi_n = NH2O*((GAMMA*GAMMA*HBAR*HBAR)/(4.*KB*Temperature));
309
         Real chi_n = NH2O*((GAMMA*GAMMA*HBAR*HBAR)/(4.*KB*Temperature));
332
         return chi_n*Porosity*B0;
310
         return chi_n*Porosity*B0;
378
                         0, step[1], step[2],
356
                         0, step[1], step[2],
379
                   step[0], step[1], step[2] ).finished();
357
                   step[0], step[1], step[2] ).finished();
380
 
358
 
381
-        //VectorXcr kvals(8);       // individual kernel vals
382
         MatrixXcr kvals(8, PulseI.size());       // individual kernel vals
359
         MatrixXcr kvals(8, PulseI.size());       // individual kernel vals
383
         cpoints->ClearFields();
360
         cpoints->ClearFields();
384
         for (int ichild=0; ichild<8; ++ichild) {
361
         for (int ichild=0; ichild<8; ++ichild) {
439
     //      Method:  EvaluateKids2 -- same as Evaluate Kids, but include VTK octree generation
416
     //      Method:  EvaluateKids2 -- same as Evaluate Kids, but include VTK octree generation
440
     //--------------------------------------------------------------------------------------
417
     //--------------------------------------------------------------------------------------
441
     void KernelV0::EvaluateKids2( const Vector3r& size, const int& level, const Vector3r& cpos,
418
     void KernelV0::EvaluateKids2( const Vector3r& size, const int& level, const Vector3r& cpos,
442
-        const Complex& parentVal, vtkHyperOctree* oct, vtkHyperOctreeCursor* curse) {
419
+        const VectorXcr& parentVal, vtkHyperOctree* oct, vtkHyperOctreeCursor* curse) {
443
 
420
 
444
         std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
421
         std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
445
         std::cout.flush();
422
         std::cout.flush();
460
                         0, step[1], step[2],
437
                         0, step[1], step[2],
461
                   step[0], step[1], step[2] ).finished();
438
                   step[0], step[1], step[2] ).finished();
462
 
439
 
463
-        VectorXcr kvals(8);                     // individual kernel vals
440
+        MatrixXcr kvals(8, PulseI.size());       // individual kernel vals
464
         cpoints->ClearFields();
441
         cpoints->ClearFields();
465
         for (int ichild=0; ichild<8; ++ichild) {
442
         for (int ichild=0; ichild<8; ++ichild) {
466
             Vector3r cp = pos;    // Eigen complains about combining these
443
             Vector3r cp = pos;    // Eigen complains about combining these
490
         }
467
         }
491
 
468
 
492
         for (int ichild=0; ichild<8; ++ichild) {
469
         for (int ichild=0; ichild<8; ++ichild) {
493
-            Vector3r cp = pos; // Eigen complains about combining these
470
+            Vector3r cp = pos;    // Eigen complains about combining these
494
             cp += posadd.row(ichild);
471
             cp += posadd.row(ichild);
495
-            kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild))(0);
472
+            kvals.row(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
496
         }
473
         }
497
 
474
 
498
-        Complex ksum = kvals.sum();     // Kernel sum
475
+        VectorXcr ksum = kvals.colwise().sum();     // Kernel sum
499
         // Evaluate whether or not furthur splitting is needed
476
         // Evaluate whether or not furthur splitting is needed
500
-        if ( std::abs(ksum - parentVal) > tol || level < minLevel && level < maxLevel ) {
501
-        //if ( std::abs(ksum.real()-parentVal.real())>tol || std::abs(ksum.imag()-parentVal.imag()) >tol ) {
477
+        if ( ((ksum - parentVal).array().abs() > tol).any() || level < minLevel && level < maxLevel ) {
502
             oct->SubdivideLeaf(curse);
478
             oct->SubdivideLeaf(curse);
503
             for (int ichild=0; ichild<8; ++ichild) {
479
             for (int ichild=0; ichild<8; ++ichild) {
504
                 curse->ToChild(ichild);
480
                 curse->ToChild(ichild);
515
                 }
491
                 }
516
                 */
492
                 */
517
                 /* End of position test */
493
                 /* End of position test */
518
-                EvaluateKids2( size, level+1, cp, kvals(ichild), oct, curse );
494
+                EvaluateKids2( size, level+1, cp, kvals.row(ichild), oct, curse );
519
                 curse->ToParent();
495
                 curse->ToParent();
520
             }
496
             }
521
             return;  // not a leaf
497
             return;  // not a leaf
522
         }
498
         }
523
-        LeafDict[curse->GetLeafId()] = ksum/(8.*vol);     //kvals(ichild) / vol;     // VTK
524
-        LeafDictIdx[curse->GetLeafId()] = nleaves;        // VTK
525
-        SUM += ksum;
499
+        LeafDict[curse->GetLeafId()] = ksum/(8.*vol);
500
+        LeafDictIdx[curse->GetLeafId()] = nleaves;
501
+        Kern.row(ilay) += ksum;
526
         VOLSUM += 8*vol;
502
         VOLSUM += 8*vol;
527
         nleaves += 1;
503
         nleaves += 1;
528
         return;     // is a leaf
504
         return;     // is a leaf

正在加载...
取消
保存