Bladeren bron

All q at once kernel evaluation for speedup, don't repeat EMEarth1D calcualtions for each moment.

master
Trevor Irons 7 jaren geleden
bovenliggende
commit
95100edf0b
3 gewijzigde bestanden met toevoegingen van 81 en 58 verwijderingen
  1. 4
    3
      examples/KernelV0.cpp
  2. 10
    9
      include/KernelV0.h
  3. 67
    46
      src/KernelV0.cpp

+ 4
- 3
examples/KernelV0.cpp Bestand weergeven

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-13 );
48
+        Kern->SetTolerance( 1e-12 );
49
         //Kern->SetTolerance( .55 ) ; // 1%
49
         //Kern->SetTolerance( .55 ) ; // 1%
50
 
50
 
51
         Kern->SetPulseDuration(0.020);
51
         Kern->SetPulseDuration(0.020);
57
              31.347205894128976, 28.06346770557137, 25.139117042424758, 22.53420773366429, 20.214205433283347,
57
              31.347205894128976, 28.06346770557137, 25.139117042424758, 22.53420773366429, 20.214205433283347,
58
              18.144318026099942, 16.299965972298878, 14.652633628829891, 13.184271405688083, 11.870540177313893,
58
              18.144318026099942, 16.299965972298878, 14.652633628829891, 13.184271405688083, 11.870540177313893,
59
              10.697057141915716, 9.64778948429609, 8.709338689612677, 7.871268012862094;
59
              10.697057141915716, 9.64778948429609, 8.709338689612677, 7.871268012862094;
60
-        Kern->SetPulseCurrent( VectorXr::LinSpaced( 1, 10, 200 )  ); // nbins, low, high
61
-        Kern->SetDepthLayerInterfaces( VectorXr::LinSpaced( 2, 5, 5.5 ) );
60
+        //Kern->SetPulseCurrent( VectorXr::LinSpaced( 1, 10, 200 )  ); // nbins, low, high
61
+        Kern->SetPulseCurrent( I ); // nbins, low, high
62
+        Kern->SetDepthLayerInterfaces( VectorXr::LinSpaced( 30, 3, 45.5 ) );
62
 
63
 
63
     // 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
64
     // may be more natural to work with?
65
     // may be more natural to work with?

+ 10
- 9
include/KernelV0.h Bestand weergeven

223
         /**
223
         /**
224
          *  Returns the kernel value for an input prism
224
          *  Returns the kernel value for an input prism
225
          */
225
          */
226
-        Complex 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);
228
+//         Complex ComputeV0Cell(const Vector3cr& Bt, const Vector3cr& Br, const Real& vol,
229
+//             const Real& phi);
230
 
230
 
231
         Complex ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
231
         Complex ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
232
                 const Real& sintheta, const Real& phase, const Real& Mn0Abs,
232
                 const Real& sintheta, const Real& phase, const Real& Mn0Abs,
236
 
236
 
237
         Vector3r ComputeMn0(const Real& Porosity, const Vector3r& B0);
237
         Vector3r ComputeMn0(const Real& Porosity, const Vector3r& B0);
238
 
238
 
239
-        Complex IntegrateOnOctreeGrid( const int& ilay, const int& iq, bool vtkOutput=false );
239
+        Complex IntegrateOnOctreeGrid( const int& iq, bool vtkOutput=false );
240
 
240
 
241
         /**
241
         /**
242
          *  Recursive call to integrate a function on an adaptive Octree Grid.
242
          *  Recursive call to integrate a function on an adaptive Octree Grid.
249
          *  @param[in] cpos is the centre position of the parent cuboid
249
          *  @param[in] cpos is the centre position of the parent cuboid
250
          */
250
          */
251
         void EvaluateKids(  const Vector3r& size, const int& level, const Vector3r& cpos,
251
         void EvaluateKids(  const Vector3r& size, const int& level, const Vector3r& cpos,
252
-                            const Complex& parentVal );
252
+                            const VectorXcr& parentVal );
253
 
253
 
254
         #ifdef LEMMAUSEVTK
254
         #ifdef LEMMAUSEVTK
255
         /**
255
         /**
264
 
264
 
265
         // ====================  DATA MEMBERS  =========================
265
         // ====================  DATA MEMBERS  =========================
266
 
266
 
267
+        int                                       ilay;
267
         int                                       nleaves;
268
         int                                       nleaves;
268
         int                                       minLevel=4;
269
         int                                       minLevel=4;
269
         int                                       maxLevel=8;
270
         int                                       maxLevel=8;
272
         Real                                      tol=1e-11;
273
         Real                                      tol=1e-11;
273
         Real                                      Temperature=283.;
274
         Real                                      Temperature=283.;
274
         Real                                      Taup = .020;  // Sec
275
         Real                                      Taup = .020;  // Sec
275
-        Real                                      Ip = 10;      // Amps
276
+        Real                                      Ip;           // Amps pulse current, deprecated PulseI
276
         Real                                      Larmor;
277
         Real                                      Larmor;
277
 
278
 
278
-        Complex                                   SUM;
279
+        Complex                                   SUM;          // Depreceated, use Kern instead
279
 
280
 
280
         Vector3r                                  Size;
281
         Vector3r                                  Size;
281
         Vector3r                                  Origin;
282
         Vector3r                                  Origin;
283
         VectorXr   PulseI;
284
         VectorXr   PulseI;
284
         VectorXr   Interfaces;
285
         VectorXr   Interfaces;
285
 
286
 
286
-        std::shared_ptr< LayeredEarthEM >         SigmaModel = nullptr;
287
+        MatrixXcr   Kern;
287
 
288
 
289
+        std::shared_ptr< LayeredEarthEM >         SigmaModel = nullptr;
288
         std::shared_ptr< FieldPoints >            cpoints;
290
         std::shared_ptr< FieldPoints >            cpoints;
289
 
291
 
290
         std::map< std::string , std::shared_ptr< PolygonalWireAntenna > >  TxRx;
292
         std::map< std::string , std::shared_ptr< PolygonalWireAntenna > >  TxRx;
291
-
292
         std::map< std::string , std::shared_ptr< EMEarth1D > >             EMEarths;
293
         std::map< std::string , std::shared_ptr< EMEarth1D > >             EMEarths;
293
 
294
 
294
         #ifdef LEMMAUSEVTK
295
         #ifdef LEMMAUSEVTK

+ 67
- 46
src/KernelV0.cpp Bestand weergeven

121
                 // TODO query for method, altough with flat antennae, this is fastest
121
                 // TODO query for method, altough with flat antennae, this is fastest
122
                 EMEarths[tx]->SetHankelTransformMethod(ANDERSON801);
122
                 EMEarths[tx]->SetHankelTransformMethod(ANDERSON801);
123
                 EMEarths[tx]->SetTxRxMode(TX);
123
                 EMEarths[tx]->SetTxRxMode(TX);
124
+                TxRx[tx]->SetCurrent(1.);
124
         }
125
         }
125
         for (auto rx : Rx) {
126
         for (auto rx : Rx) {
126
             if (EMEarths.count(rx)) {
127
             if (EMEarths.count(rx)) {
134
                     // TODO query for method, altough with flat antennae, this is fastest
135
                     // TODO query for method, altough with flat antennae, this is fastest
135
                     EMEarths[rx]->SetHankelTransformMethod(ANDERSON801);
136
                     EMEarths[rx]->SetHankelTransformMethod(ANDERSON801);
136
                     EMEarths[rx]->SetTxRxMode(RX);
137
                     EMEarths[rx]->SetTxRxMode(RX);
138
+                    TxRx[rx]->SetCurrent(1.);
137
             }
139
             }
138
         }
140
         }
139
 
141
 
140
         std::cout << "Calculating K0 kernel\n";
142
         std::cout << "Calculating K0 kernel\n";
141
-        MatrixXcr Kern = MatrixXcr::Zero( Interfaces.size()-1, PulseI.size() );
142
-        for (int ilay=0; ilay<Interfaces.size()-1; ++ilay) {
143
-            for (int iq=0; iq< PulseI.size(); ++iq) {
144
-                std::cout << "Layer " << ilay << " q " << iq << std::endl;
145
-                Size(2) = Interfaces(ilay+1) - Interfaces(ilay);
146
-                Origin(2) = Interfaces(ilay);
147
-                Ip = PulseI(iq);
148
-                Kern(ilay, iq) = IntegrateOnOctreeGrid( ilay, iq, vtkOutput );
149
-            }
143
+        Kern = MatrixXcr::Zero( Interfaces.size()-1, PulseI.size() );
144
+        for (ilay=0; ilay<Interfaces.size()-1; ++ilay) {
145
+            std::cout << "Layer " << ilay << std::endl; //<< " q " << iq << std::endl;
146
+            Size(2) = Interfaces(ilay+1) - 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
+            //}
150
         }
153
         }
151
         std::cout << "\rFinished KERNEL\n";
154
         std::cout << "\rFinished KERNEL\n";
152
         std::cout << "real\n";
155
         std::cout << "real\n";
153
         std::cout << Kern.real() << std::endl;
156
         std::cout << Kern.real() << std::endl;
154
         std::cout << "imag\n";
157
         std::cout << "imag\n";
155
         std::cout << Kern.imag() << std::endl;
158
         std::cout << Kern.imag() << std::endl;
156
-        //IntegrateOnOctreeGrid( vtkOutput );
157
     }
159
     }
158
 
160
 
159
     //--------------------------------------------------------------------------------------
161
     //--------------------------------------------------------------------------------------
160
     //       Class:  KernelV0
162
     //       Class:  KernelV0
161
     //      Method:  IntegrateOnOctreeGrid
163
     //      Method:  IntegrateOnOctreeGrid
162
     //--------------------------------------------------------------------------------------
164
     //--------------------------------------------------------------------------------------
163
-    Complex KernelV0::IntegrateOnOctreeGrid( const int& ilay, const int& iq, bool vtkOutput) {
165
+    Complex KernelV0::IntegrateOnOctreeGrid( const int& iq, bool vtkOutput) {
164
 
166
 
165
         Vector3r cpos = Origin + Size/2.;
167
         Vector3r cpos = Origin + Size/2.;
166
 
168
 
168
         VOLSUM = 0;
170
         VOLSUM = 0;
169
         nleaves = 0;
171
         nleaves = 0;
170
         if (!vtkOutput) {
172
         if (!vtkOutput) {
171
-            EvaluateKids( Size, 0, cpos, 1e6 );
173
+            EvaluateKids( Size, 0, cpos, VectorXcr::Ones(PulseI.size()) );
172
         } else {
174
         } else {
173
         #ifdef LEMMAUSEVTK
175
         #ifdef LEMMAUSEVTK
174
             vtkHyperOctree* oct = vtkHyperOctree::New();
176
             vtkHyperOctree* oct = vtkHyperOctree::New();
253
     //       Class:  KernelV0
255
     //       Class:  KernelV0
254
     //      Method:  f
256
     //      Method:  f
255
     //--------------------------------------------------------------------------------------
257
     //--------------------------------------------------------------------------------------
256
-    Complex KernelV0::f( const Vector3r& r, const Real& volume, const Vector3cr& Ht, const Vector3cr& Hr ) {
257
-        //return Complex(volume*Ht.dot(Hr));
258
-        return ComputeV0Cell(MU0*Ht, MU0*Hr, volume, 1.0);
259
-    }
260
-
261
-    //--------------------------------------------------------------------------------------
262
-    //       Class:  KernelV0
263
-    //      Method:  ComputeV0Cell
264
-    //--------------------------------------------------------------------------------------
265
-    Complex KernelV0::ComputeV0Cell(const Vector3cr& Bt,
266
-                const Vector3cr& Br, const Real& vol, const Real& phi) {
258
+    VectorXcr KernelV0::f( const Vector3r& r, const Real& volume, const Vector3cr& Ht, const Vector3cr& Hr ) {
267
 
259
 
268
         // Compute the elliptic fields
260
         // Compute the elliptic fields
269
         Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
261
         Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
270
         Vector3r B0 = SigmaModel->GetMagneticField();
262
         Vector3r B0 = SigmaModel->GetMagneticField();
271
 
263
 
272
         // Elliptic representation
264
         // Elliptic representation
273
-        EllipticB EBT = EllipticFieldRep(Bt, B0hat);
274
-        EllipticB EBR = EllipticFieldRep(Br, B0hat);
265
+        EllipticB EBT = EllipticFieldRep(MU0*Ht, B0hat);
266
+        EllipticB EBR = EllipticFieldRep(MU0*Hr, B0hat);
275
 
267
 
276
         // Compute Mn0
268
         // Compute Mn0
277
-        Vector3r Mn0 = ComputeMn0(phi, B0);
269
+        Vector3r Mn0 = ComputeMn0(1.0, B0);
278
         Real Mn0Abs = Mn0.norm();
270
         Real Mn0Abs = Mn0.norm();
279
 
271
 
280
-        // Compute the tipping angle
281
-        Real sintheta = std::sin(0.5*GAMMA*Ip*Taup*std::abs(EBT.alpha-EBT.beta));
282
-
283
-        // Compute phase delay, TODO add transmiiter phase and delay time phase!
272
+        // Compute phase delay
273
+        // TODO add transmiiter phase and delay time phase!
284
         Real phase = EBR.zeta+EBT.zeta;
274
         Real phase = EBR.zeta+EBT.zeta;
285
 
275
 
286
-        return ComputeV0Cell(EBT, EBR, sintheta, phase, Mn0Abs, vol);
276
+        // Calcuate vector of all responses
277
+        VectorXcr F = VectorXcr::Zero( PulseI.size() );
278
+        for (int iq=0; iq<PulseI.size(); ++iq) {
279
+            // Compute the tipping angle
280
+            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);
282
+        }
283
+        return F;
287
     }
284
     }
288
 
285
 
286
+//     //--------------------------------------------------------------------------------------
287
+//     //       Class:  KernelV0
288
+//     //      Method:  ComputeV0Cell
289
+//     //--------------------------------------------------------------------------------------
290
+//     Complex KernelV0::ComputeV0Cell(const Vector3cr& Bt,
291
+//                 const Vector3cr& Br, const Real& vol, const Real& phi) {
292
+//
293
+//         // Compute the elliptic fields
294
+//         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);
310
+//     }
311
+
289
     //--------------------------------------------------------------------------------------
312
     //--------------------------------------------------------------------------------------
290
     //       Class:  KernelV0
313
     //       Class:  KernelV0
291
     //      Method:  ComputeV0Cell
314
     //      Method:  ComputeV0Cell
293
     Complex KernelV0::ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
316
     Complex KernelV0::ComputeV0Cell(const EllipticB& EBT, const EllipticB& EBR,
294
                 const Real& sintheta, const Real& phase, const Real& Mn0Abs,
317
                 const Real& sintheta, const Real& phase, const Real& Mn0Abs,
295
                 const Real& vol) {
318
                 const Real& vol) {
296
-
297
-        Vector3r B0hat = {1,0,0};
298
-
299
         // earth response of receiver adjoint field
319
         // earth response of receiver adjoint field
320
+        Vector3r B0hat = SigmaModel->GetMagneticFieldUnitVector();
300
         Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
321
         Complex ejztr = std::exp(Complex(0, EBR.zeta + EBT.zeta));
301
-
302
-        Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) +
303
-               (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
322
+        Complex PhaseTerm = EBR.bhat.dot(EBT.bhat) + (B0hat.dot(EBR.bhat.cross(EBT.bhat) ));
304
         return -vol*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
323
         return -vol*Complex(0,Larmor)*Mn0Abs*(EBR.alpha+EBR.beta)*ejztr*sintheta*PhaseTerm;
305
     }
324
     }
306
 
325
 
338
     //      Method:  EvaluateKids
357
     //      Method:  EvaluateKids
339
     //--------------------------------------------------------------------------------------
358
     //--------------------------------------------------------------------------------------
340
     void KernelV0::EvaluateKids( const Vector3r& size, const int& level, const Vector3r& cpos,
359
     void KernelV0::EvaluateKids( const Vector3r& size, const int& level, const Vector3r& cpos,
341
-        const Complex& parentVal ) {
360
+        const VectorXcr& parentVal ) {
342
 
361
 
343
         std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
362
         std::cout << "\r" << (int)(1e2*VOLSUM/(Size[0]*Size[1]*Size[2])) << "\t" << nleaves;
344
         std::cout.flush();
363
         std::cout.flush();
359
                         0, step[1], step[2],
378
                         0, step[1], step[2],
360
                   step[0], step[1], step[2] ).finished();
379
                   step[0], step[1], step[2] ).finished();
361
 
380
 
362
-        VectorXcr kvals(8);       // individual kernel vals
381
+        //VectorXcr kvals(8);       // individual kernel vals
382
+        MatrixXcr kvals(8, PulseI.size());       // individual kernel vals
363
         cpoints->ClearFields();
383
         cpoints->ClearFields();
364
         for (int ichild=0; ichild<8; ++ichild) {
384
         for (int ichild=0; ichild<8; ++ichild) {
365
             Vector3r cp = pos;    // Eigen complains about combining these
385
             Vector3r cp = pos;    // Eigen complains about combining these
392
         for (int ichild=0; ichild<8; ++ichild) {
412
         for (int ichild=0; ichild<8; ++ichild) {
393
             Vector3r cp = pos;    // Eigen complains about combining these
413
             Vector3r cp = pos;    // Eigen complains about combining these
394
             cp += posadd.row(ichild);
414
             cp += posadd.row(ichild);
395
-            kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
415
+            kvals.row(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
396
         }
416
         }
397
 
417
 
398
-        Complex ksum = kvals.sum();     // Kernel sum
418
+        VectorXcr ksum = kvals.colwise().sum();     // Kernel sum
399
         // Evaluate whether or not furthur splitting is needed
419
         // Evaluate whether or not furthur splitting is needed
400
-        if ( std::abs(ksum - parentVal) > tol || level < minLevel && level < maxLevel ) {
420
+        if ( ((ksum - parentVal).array().abs() > tol).any() || level < minLevel && level < maxLevel ) {
421
+            // Not a leaf dive further in
401
             for (int ichild=0; ichild<8; ++ichild) {
422
             for (int ichild=0; ichild<8; ++ichild) {
402
                 Vector3r cp = pos; // Eigen complains about combining these
423
                 Vector3r cp = pos; // Eigen complains about combining these
403
                 cp += posadd.row(ichild);
424
                 cp += posadd.row(ichild);
404
-                EvaluateKids( size, level+1, cp, kvals(ichild) );
425
+                EvaluateKids( size, level+1, cp, kvals.row(ichild) );
405
             }
426
             }
406
             return; // not leaf
427
             return; // not leaf
407
         }
428
         }
408
-        // Save here instead?
409
-        SUM += ksum;
429
+        // implicit else, is a leaf
430
+        Kern.row(ilay) += ksum;
410
         VOLSUM += 8.*vol;
431
         VOLSUM += 8.*vol;
411
         nleaves += 1;
432
         nleaves += 1;
412
         return;     // is leaf
433
         return;     // is leaf
471
         for (int ichild=0; ichild<8; ++ichild) {
492
         for (int ichild=0; ichild<8; ++ichild) {
472
             Vector3r cp = pos; // Eigen complains about combining these
493
             Vector3r cp = pos; // Eigen complains about combining these
473
             cp += posadd.row(ichild);
494
             cp += posadd.row(ichild);
474
-            kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild));
495
+            kvals(ichild) = f(cp, vol, Ht.col(ichild), Hr.col(ichild))(0);
475
         }
496
         }
476
 
497
 
477
         Complex ksum = kvals.sum();     // Kernel sum
498
         Complex ksum = kvals.sum();     // Kernel sum

Laden…
Annuleren
Opslaan