GEOS
PrecomputeSourcesAndReceiversKernel.hpp
Go to the documentation of this file.
1 /*
2  * ------------------------------------------------------------------------------------------------------------
3  * SPDX-License-Identifier: LGPL-2.1-only
4  *
5  * Copyright (c) 2016-2024 Lawrence Livermore National Security LLC
6  * Copyright (c) 2018-2024 TotalEnergies
7  * Copyright (c) 2018-2024 The Board of Trustees of the Leland Stanford Junior University
8  * Copyright (c) 2023-2024 Chevron
9  * Copyright (c) 2019- GEOS/GEOSX Contributors
10  * All rights reserved
11  *
12  * See top level LICENSE, COPYRIGHT, CONTRIBUTORS, NOTICE, and ACKNOWLEDGEMENTS files for details.
13  * ------------------------------------------------------------------------------------------------------------
14  */
15 
20 #ifndef GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_PRECOMPUTESOURCESANDRECEIVERSKERNEL_HPP_
21 #define GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_PRECOMPUTESOURCESANDRECEIVERSKERNEL_HPP_
22 
23 namespace geos
24 {
25 
27 {
28 
29  using EXEC_POLICY = parallelDevicePolicy< >;
30 
55  template< typename EXEC_POLICY, typename FE_TYPE >
56  static void
58  ArrayOfArraysView< localIndex const > const baseFacesToNodes,
60  arrayView1d< globalIndex const > const baseNodeLocalToGlobal,
61  arrayView1d< globalIndex const > const elementLocalToGlobal,
62  ArrayOfArraysView< localIndex const > const baseNodesToElements,
64  arrayView1d< integer const > const elemGhostRank,
66  arrayView2d< localIndex const > const elemsToFaces,
67  arrayView2d< real64 const > const & elemCenter,
68  arrayView2d< real64 const > const sourceCoordinates,
69  arrayView1d< localIndex > const sourceIsAccessible,
70  arrayView2d< localIndex > const sourceNodeIds,
71  arrayView2d< real64 > const sourceConstants,
72  arrayView2d< real64 const > const receiverCoordinates,
73  arrayView1d< localIndex > const receiverIsLocal,
74  arrayView2d< localIndex > const receiverNodeIds,
75  arrayView2d< real64 > const receiverConstants,
76  arrayView2d< real32 > const sourceValue,
77  real64 const dt,
78  real32 const timeSourceFrequency,
79  real32 const timeSourceDelay,
80  localIndex const rickerOrder,
81  arrayView1d< TableFunction::KernelWrapper const > const sourceWaveletTableWrappers,
82  bool useSourceWaveletTables )
83  {
84  constexpr localIndex numNodesPerElem = FE_TYPE::numNodes;
85 
86  forAll< EXEC_POLICY >( size, [=] GEOS_HOST_DEVICE ( localIndex const k )
87  {
88  real64 const center[3] = { elemCenter[k][0],
89  elemCenter[k][1],
90  elemCenter[k][2] };
91 
92  // Step 1: locate the sources, and precompute the source term
93 
95  for( localIndex isrc = 0; isrc < sourceCoordinates.size( 0 ); ++isrc )
96  {
97  if( sourceIsAccessible[isrc] == 0 )
98  {
99  real64 const coords[3] = { sourceCoordinates[isrc][0],
100  sourceCoordinates[isrc][1],
101  sourceCoordinates[isrc][2] };
102  bool const sourceFound =
103  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
104  baseNodeCoords,
105  elemsToFaces,
106  baseFacesToNodes,
107  baseNodesToElements,
108  baseNodeLocalToGlobal,
109  elementLocalToGlobal,
110  center,
111  coords );
112  if( sourceFound )
113  {
114  real64 coordsOnRefElem[3]{};
115 
116 
117  WaveSolverUtils::computeCoordinatesOnReferenceElement< FE_TYPE >( coords,
118  baseElemsToNodes[k],
119  baseNodeCoords,
120  coordsOnRefElem );
121 
122  sourceIsAccessible[isrc] = 1;
123  real64 Ntest[numNodesPerElem];
124  FE_TYPE::calcN( coordsOnRefElem, Ntest );
125 
126  for( localIndex a = 0; a < numNodesPerElem; ++a )
127  {
128  sourceNodeIds[isrc][a] = elemsToNodes( k, a );
129  sourceConstants[isrc][a] = Ntest[a];
130  }
131 
132  for( localIndex cycle = 0; cycle < sourceValue.size( 0 ); ++cycle )
133  {
134  real64 const time_n = cycle * dt;
135  if( useSourceWaveletTables )
136  {
137  sourceValue[cycle][isrc]= sourceWaveletTableWrappers[ isrc ].compute( &time_n );
138  }
139  else
140  {
141  sourceValue[cycle][isrc] = WaveSolverUtils::evaluateRicker( cycle * dt, timeSourceFrequency, timeSourceDelay, rickerOrder );
142  }
143  }
144 
145  }
146  }
147  } // end loop over all sources
148 
149 
150  // Step 2: locate the receivers, and precompute the receiver term
151 
153  for( localIndex ircv = 0; ircv < receiverCoordinates.size( 0 ); ++ircv )
154  {
155  if( receiverIsLocal[ircv] == 0 )
156  {
157  real64 const coords[3] = { receiverCoordinates[ircv][0],
158  receiverCoordinates[ircv][1],
159  receiverCoordinates[ircv][2] };
160 
161  real64 coordsOnRefElem[3]{};
162  bool const receiverFound =
163  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
164  baseNodeCoords,
165  elemsToFaces,
166  baseFacesToNodes,
167  baseNodesToElements,
168  baseNodeLocalToGlobal,
169  elementLocalToGlobal,
170  center,
171  coords );
172 
173  if( receiverFound && elemGhostRank[k] < 0 )
174  {
175  WaveSolverUtils::computeCoordinatesOnReferenceElement< FE_TYPE >( coords,
176  baseElemsToNodes[k],
177  baseNodeCoords,
178  coordsOnRefElem );
179 
180  receiverIsLocal[ircv] = 1;
181 
182  real64 Ntest[numNodesPerElem];
183  FE_TYPE::calcN( coordsOnRefElem, Ntest );
184 
185  for( localIndex a = 0; a < numNodesPerElem; ++a )
186  {
187  receiverNodeIds[ircv][a] = elemsToNodes( k, a );
188  receiverConstants[ircv][a] = Ntest[a];
189  }
190  }
191  }
192  } // end loop over receivers
193 
194  } );
195 
196  }
197 
198 
226  template< typename EXEC_POLICY, typename FE_TYPE >
227  static void
229  localIndex const regionIndex,
230  ArrayOfArraysView< localIndex const > const baseFacesToNodes,
232  arrayView1d< globalIndex const > const baseNodeLocalToGlobal,
233  arrayView1d< globalIndex const > const elementLocalToGlobal,
234  ArrayOfArraysView< localIndex const > const baseNodesToElements,
236  arrayView1d< integer const > const elemGhostRank,
238  arrayView2d< localIndex const > const elemsToFaces,
239  arrayView2d< real64 const > const & elemCenter,
240  arrayView2d< real64 const > const sourceCoordinates,
241  arrayView1d< localIndex > const sourceIsAccessible,
242  arrayView1d< localIndex > const sourceElem,
243  arrayView2d< localIndex > const sourceNodeIds,
244  arrayView2d< real64 > const sourceConstants,
245  arrayView1d< localIndex > const sourceRegion,
246  arrayView2d< real64 const > const receiverCoordinates,
247  arrayView1d< localIndex > const receiverIsLocal,
248  arrayView1d< localIndex > const receiverElem,
249  arrayView2d< localIndex > const receiverNodeIds,
250  arrayView2d< real64 > const receiverConstants,
251  arrayView1d< localIndex > const receiverRegion )
252  {
253  constexpr localIndex numNodesPerElem = FE_TYPE::numNodes;
254 
255  forAll< EXEC_POLICY >( size, [=] GEOS_HOST_DEVICE ( localIndex const k )
256  {
257  real64 const center[3] = { elemCenter[k][0],
258  elemCenter[k][1],
259  elemCenter[k][2] };
260 
261  // Step 1: locate the sources, and precompute the source term
262 
264  for( localIndex isrc = 0; isrc < sourceCoordinates.size( 0 ); ++isrc )
265  {
266  if( sourceIsAccessible[isrc] == 0 )
267  {
268  real64 const coords[3] = { sourceCoordinates[isrc][0],
269  sourceCoordinates[isrc][1],
270  sourceCoordinates[isrc][2] };
271 
272  bool const sourceFound =
273  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
274  baseNodeCoords,
275  elemsToFaces,
276  baseFacesToNodes,
277  baseNodesToElements,
278  baseNodeLocalToGlobal,
279  elementLocalToGlobal,
280  center,
281  coords );
282  if( sourceFound )
283  {
284  real64 coordsOnRefElem[3]{};
285 
286  WaveSolverUtils::computeCoordinatesOnReferenceElement< FE_TYPE >( coords,
287  baseElemsToNodes[k],
288  baseNodeCoords,
289  coordsOnRefElem );
290 
291  sourceIsAccessible[isrc] = 1;
292  sourceElem[isrc] = k;
293  sourceRegion[isrc] = regionIndex;
294  real64 Ntest[numNodesPerElem];
295  FE_TYPE::calcN( coordsOnRefElem, Ntest );
296 
297  for( localIndex a = 0; a < numNodesPerElem; ++a )
298  {
299  sourceNodeIds[isrc][a] = elemsToNodes[k][a];
300  sourceConstants[isrc][a] = Ntest[a];
301  }
302 
303  }
304  }
305  } // end loop over all sources
306 
307 
308  // Step 2: locate the receivers, and precompute the receiver term
309 
311  for( localIndex ircv = 0; ircv < receiverCoordinates.size( 0 ); ++ircv )
312  {
313  if( receiverIsLocal[ircv] == 0 )
314  {
315  real64 const coords[3] = { receiverCoordinates[ircv][0],
316  receiverCoordinates[ircv][1],
317  receiverCoordinates[ircv][2] };
318 
319  real64 coordsOnRefElem[3]{};
320  bool const receiverFound =
321  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
322  baseNodeCoords,
323  elemsToFaces,
324  baseFacesToNodes,
325  baseNodesToElements,
326  baseNodeLocalToGlobal,
327  elementLocalToGlobal,
328  center,
329  coords );
330 
331  if( receiverFound && elemGhostRank[k] < 0 )
332  {
333  WaveSolverUtils::computeCoordinatesOnReferenceElement< FE_TYPE >( coords,
334  baseElemsToNodes[k],
335  baseNodeCoords,
336  coordsOnRefElem );
337  receiverIsLocal[ircv] = 1;
338  receiverElem[ircv] = k;
339  receiverRegion[ircv] = regionIndex;
340 
341  real64 Ntest[numNodesPerElem];
342  FE_TYPE::calcN( coordsOnRefElem, Ntest );
343 
344  for( localIndex a = 0; a < numNodesPerElem; ++a )
345  {
346  receiverNodeIds[ircv][a] = elemsToNodes[k][a];
347  receiverConstants[ircv][a] = Ntest[a];
348  }
349  }
350  }
351  } // end loop over receivers
352 
353  } );
354 
355  }
356 
357 
358 
392  template< typename EXEC_POLICY, typename FE_TYPE >
393  static void
395  ArrayOfArraysView< localIndex const > const baseFacesToNodes,
397  arrayView1d< globalIndex const > const baseNodeLocalToGlobal,
398  arrayView1d< globalIndex const > const elementLocalToGlobal,
399  ArrayOfArraysView< localIndex const > const baseNodesToElements,
401  arrayView1d< integer const > const elemGhostRank,
403  arrayView2d< localIndex const > const elemsToFaces,
404  arrayView2d< real64 const > const & elemCenter,
405  arrayView2d< real64 const > const sourceCoordinates,
406  arrayView1d< localIndex > const sourceIsAccessible,
407  arrayView2d< localIndex > const sourceNodeIds,
408  arrayView2d< real64 > const sourceConstantsx,
409  arrayView2d< real64 > const sourceConstantsy,
410  arrayView2d< real64 > const sourceConstantsz,
411  arrayView2d< real64 const > const receiverCoordinates,
412  arrayView1d< localIndex > const receiverIsLocal,
413  arrayView2d< localIndex > const receiverNodeIds,
414  arrayView2d< real64 > const receiverConstants,
416  integer linearDASSamples,
417  arrayView2d< real64 const > const linearDASGeometry,
418  R1Tensor const sourceForce,
419  R2SymTensor const sourceMoment )
420  {
421  constexpr localIndex numNodesPerElem = FE_TYPE::numNodes;
422  integer nSamples = useDAS == WaveSolverUtils::DASType::none ? 1 : linearDASSamples;
423  array1d< real64 > const samplePointLocationsA( nSamples );
424  arrayView1d< real64 > const samplePointLocations = samplePointLocationsA.toView();
425  array1d< real64 > const sampleIntegrationConstantsA( nSamples );
426  arrayView1d< real64 > const sampleIntegrationConstants = sampleIntegrationConstantsA.toView();
427 
428  forAll< EXEC_POLICY >( size, [=] GEOS_HOST_DEVICE ( localIndex const k )
429  {
430  real64 const center[3] = { elemCenter[k][0],
431  elemCenter[k][1],
432  elemCenter[k][2] };
433 
434  // Step 1: locate the sources, and precompute the source term
435 
437  for( localIndex isrc = 0; isrc < sourceCoordinates.size( 0 ); ++isrc )
438  {
439  if( sourceIsAccessible[isrc] == 0 )
440  {
441  real64 const coords[3] = { sourceCoordinates[isrc][0],
442  sourceCoordinates[isrc][1],
443  sourceCoordinates[isrc][2] };
444 
445  real64 xLocal[8][3];
446 
447  for( localIndex a = 0; a < 8; ++a )
448  {
449  for( localIndex i = 0; i < 3; ++i )
450  {
451  xLocal[a][i] = baseNodeCoords( baseElemsToNodes( k, a ), i );
452  }
453  }
454 
455 
456  bool const sourceFound =
457  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
458  baseNodeCoords,
459  elemsToFaces,
460  baseFacesToNodes,
461  baseNodesToElements,
462  baseNodeLocalToGlobal,
463  elementLocalToGlobal,
464  center,
465  coords );
466 
467  if( sourceFound )
468  {
469  real64 coordsOnRefElem[3]{};
470 
471 
472  WaveSolverUtils::computeCoordinatesOnReferenceElement< FE_TYPE >( coords,
473  baseElemsToNodes[k],
474  baseNodeCoords,
475  coordsOnRefElem );
476  sourceIsAccessible[isrc] = 1;
477 
478  real64 N[numNodesPerElem];
479  real64 gradN[numNodesPerElem][3];
480  FE_TYPE::calcN( coordsOnRefElem, N );
481  FE_TYPE::calcGradNWithCorners( coordsOnRefElem, xLocal, gradN );
482  R2SymTensor moment = sourceMoment;
483  for( localIndex q=0; q< numNodesPerElem; ++q )
484  {
485  real64 inc[3] = { 0, 0, 0 };
486  sourceNodeIds[isrc][q] = elemsToNodes( k, q );
487  inc[0] += sourceForce[0] * N[q];
488  inc[1] += sourceForce[1] * N[q];
489  inc[2] += sourceForce[2] * N[q];
490 
491  LvArray::tensorOps::Ri_add_symAijBj< 3 >( inc, moment.data, gradN[q] );
492  sourceConstantsx[isrc][q] += inc[0];
493  sourceConstantsy[isrc][q] += inc[1];
494  sourceConstantsz[isrc][q] += inc[2];
495  }
496 
497  }
498  }
499  } // end loop over all sources
500 
501  // Step 2: locate the receivers, and precompute the receiver term
502 
503  // for geophones, we need only a point per receiver.
504  // for DAS, we need multiple points
505 
507  if( nSamples == 1 )
508  {
509  samplePointLocations[ 0 ] = 0;
510  }
511  else
512  {
513  for( integer i = 0; i < nSamples; ++i )
514  {
515  samplePointLocations[ i ] = -0.5 + (real64) i / ( linearDASSamples - 1 );
516  }
517  }
518 
521  if( useDAS == WaveSolverUtils::DASType::dipole )
522  {
523  sampleIntegrationConstants[ 0 ] = -1.0;
524  sampleIntegrationConstants[ 1 ] = 1.0;
525  }
527  else if( nSamples == 1 )
528  {
529  sampleIntegrationConstants[ 0 ] = 1.0;
530  }
531  else
532  {
533  for( integer i = 0; i < linearDASSamples; i++ )
534  {
535  sampleIntegrationConstants[ i ] = 1.0 / nSamples;
536  }
537  }
538 
540  for( localIndex ircv = 0; ircv < receiverCoordinates.size( 0 ); ++ircv )
541  {
542  R1Tensor receiverCenter = { receiverCoordinates[ ircv ][ 0 ], receiverCoordinates[ ircv ][ 1 ], receiverCoordinates[ ircv ][ 2 ] };
543  R1Tensor receiverVector;
544  if( useDAS == WaveSolverUtils::DASType::none )
545  {
546  receiverVector = { 0, 0, 0 };
547  }
548  else
549  {
550  receiverVector = WaveSolverUtils::computeDASVector( linearDASGeometry[ ircv ][ 0 ], linearDASGeometry[ ircv ][ 1 ] );
551  }
552  real64 receiverLength = useDAS == WaveSolverUtils::DASType::none ? 0 : linearDASGeometry[ ircv ][ 2 ];
554  for( integer iSample = 0; iSample < nSamples; ++iSample )
555  {
557  real64 const coords[3] = { receiverCenter[ 0 ] + receiverVector[ 0 ] * receiverLength * samplePointLocations[ iSample ],
558  receiverCenter[ 1 ] + receiverVector[ 1 ] * receiverLength * samplePointLocations[ iSample ],
559  receiverCenter[ 2 ] + receiverVector[ 2 ] * receiverLength * samplePointLocations[ iSample ] };
560  bool const sampleFound =
561  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
562  baseNodeCoords,
563  elemsToFaces,
564  baseFacesToNodes,
565  baseNodesToElements,
566  baseNodeLocalToGlobal,
567  elementLocalToGlobal,
568  center,
569  coords );
570  if( sampleFound && elemGhostRank[k] < 0 )
571  {
572  real64 coordsOnRefElem[3]{};
573  real64 xLocal[8][3];
574 
575  for( localIndex a = 0; a < 8; ++a )
576  {
577  for( localIndex i=0; i < 3; ++i )
578  {
579  xLocal[a][i] = baseNodeCoords( baseElemsToNodes( k, a ), i );
580  }
581  }
582 
583  WaveSolverUtils::computeCoordinatesOnReferenceElement< FE_TYPE >( coords,
584  baseElemsToNodes[k],
585  baseNodeCoords,
586  coordsOnRefElem );
587  real64 N[numNodesPerElem];
588  real64 gradN[numNodesPerElem][3];
589  FE_TYPE::calcN( coordsOnRefElem, N );
590  FE_TYPE::calcGradNWithCorners( coordsOnRefElem, xLocal, gradN );
591  for( localIndex a = 0; a < numNodesPerElem; ++a )
592  {
593  receiverNodeIds[ircv][iSample * numNodesPerElem + a] = elemsToNodes( k,
594  a );
596  {
597  receiverConstants[ircv][iSample * numNodesPerElem + a] += ( gradN[a][0] * receiverVector[0] + gradN[a][1] * receiverVector[1] + gradN[a][2] * receiverVector[2] ) *
598  sampleIntegrationConstants[ iSample ];
599  }
600  else
601  {
602  receiverConstants[ircv][iSample * numNodesPerElem + a] += N[a] * sampleIntegrationConstants[ iSample ];
603  }
604  }
605  receiverIsLocal[ ircv ] = 2;
606  }
607  } // end loop over samples
608  // determine if the current rank is the owner of this receiver
609  real64 const coords[3] = { receiverCenter[ 0 ], receiverCenter[ 1 ], receiverCenter[ 2 ] };
610  bool const receiverFound =
611  computationalGeometry::isPointInsideConvexPolyhedronRobust( k,
612  baseNodeCoords,
613  elemsToFaces,
614  baseFacesToNodes,
615  baseNodesToElements,
616  baseNodeLocalToGlobal,
617  elementLocalToGlobal,
618  center,
619  coords );
620  if( receiverFound && elemGhostRank[k] < 0 )
621  {
622  receiverIsLocal[ ircv ] = 1;
623  }
624  } // end loop over receivers
625  } );
626 
627  }
628 
629 };
630 
631 } // namespace geos
632 
633 #endif //GEOS_PHYSICSSOLVERS_WAVEPROPAGATION_PRECOMPUTESOURCESANDRECEIVERSKERNEL_HPP_
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
T data[SIZE]
Underlying array.
Definition: Tensor.hpp:142
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:180
float real32
32-bit floating point type.
Definition: DataTypes.hpp:97
LvArray::ArrayOfArraysView< T, INDEX_TYPE const, CONST_SIZES, LvArray::ChaiBuffer > ArrayOfArraysView
View of array of variable-sized arrays. See LvArray::ArrayOfArraysView for details.
Definition: DataTypes.hpp:286
double real64
64-bit floating point type.
Definition: DataTypes.hpp:99
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:85
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:196
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:176
static void Compute1DSourceAndReceiverConstants(localIndex const size, ArrayOfArraysView< localIndex const > const baseFacesToNodes, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const baseNodeCoords, arrayView1d< globalIndex const > const baseNodeLocalToGlobal, arrayView1d< globalIndex const > const elementLocalToGlobal, ArrayOfArraysView< localIndex const > const baseNodesToElements, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &baseElemsToNodes, arrayView1d< integer const > const elemGhostRank, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &elemsToNodes, arrayView2d< localIndex const > const elemsToFaces, arrayView2d< real64 const > const &elemCenter, arrayView2d< real64 const > const sourceCoordinates, arrayView1d< localIndex > const sourceIsAccessible, arrayView2d< localIndex > const sourceNodeIds, arrayView2d< real64 > const sourceConstants, arrayView2d< real64 const > const receiverCoordinates, arrayView1d< localIndex > const receiverIsLocal, arrayView2d< localIndex > const receiverNodeIds, arrayView2d< real64 > const receiverConstants, arrayView2d< real32 > const sourceValue, real64 const dt, real32 const timeSourceFrequency, real32 const timeSourceDelay, localIndex const rickerOrder, arrayView1d< TableFunction::KernelWrapper const > const sourceWaveletTableWrappers, bool useSourceWaveletTables)
Launches the precomputation of the source and receiver terms for 1D solution (2nd order acoustic)
static void Compute1DSourceAndReceiverConstantsWithElementsAndRegionStorage(localIndex const size, localIndex const regionIndex, ArrayOfArraysView< localIndex const > const baseFacesToNodes, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const baseNodeCoords, arrayView1d< globalIndex const > const baseNodeLocalToGlobal, arrayView1d< globalIndex const > const elementLocalToGlobal, ArrayOfArraysView< localIndex const > const baseNodesToElements, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &baseElemsToNodes, arrayView1d< integer const > const elemGhostRank, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &elemsToNodes, arrayView2d< localIndex const > const elemsToFaces, arrayView2d< real64 const > const &elemCenter, arrayView2d< real64 const > const sourceCoordinates, arrayView1d< localIndex > const sourceIsAccessible, arrayView1d< localIndex > const sourceElem, arrayView2d< localIndex > const sourceNodeIds, arrayView2d< real64 > const sourceConstants, arrayView1d< localIndex > const sourceRegion, arrayView2d< real64 const > const receiverCoordinates, arrayView1d< localIndex > const receiverIsLocal, arrayView1d< localIndex > const receiverElem, arrayView2d< localIndex > const receiverNodeIds, arrayView2d< real64 > const receiverConstants, arrayView1d< localIndex > const receiverRegion)
Launches the precomputation of the source and receiver terms with storage of elements and region in w...
static void Compute3DSourceAndReceiverConstantsWithDAS(localIndex const size, ArrayOfArraysView< localIndex const > const baseFacesToNodes, arrayView2d< real64 const, nodes::REFERENCE_POSITION_USD > const baseNodeCoords, arrayView1d< globalIndex const > const baseNodeLocalToGlobal, arrayView1d< globalIndex const > const elementLocalToGlobal, ArrayOfArraysView< localIndex const > const baseNodesToElements, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &baseElemsToNodes, arrayView1d< integer const > const elemGhostRank, arrayView2d< localIndex const, cells::NODE_MAP_USD > const &elemsToNodes, arrayView2d< localIndex const > const elemsToFaces, arrayView2d< real64 const > const &elemCenter, arrayView2d< real64 const > const sourceCoordinates, arrayView1d< localIndex > const sourceIsAccessible, arrayView2d< localIndex > const sourceNodeIds, arrayView2d< real64 > const sourceConstantsx, arrayView2d< real64 > const sourceConstantsy, arrayView2d< real64 > const sourceConstantsz, arrayView2d< real64 const > const receiverCoordinates, arrayView1d< localIndex > const receiverIsLocal, arrayView2d< localIndex > const receiverNodeIds, arrayView2d< real64 > const receiverConstants, WaveSolverUtils::DASType useDAS, integer linearDASSamples, arrayView2d< real64 const > const linearDASGeometry, R1Tensor const sourceForce, R2SymTensor const sourceMoment)
Launches the precomputation of the source and receiver terms for 3D arrays solution and DAS receiver ...
static GEOS_HOST_DEVICE R1Tensor computeDASVector(real64 const dip, real64 const azimuth)
Converts the DAS direction from dip/azimuth to a 3D unit vector.
@ none
deactivate DAS computation
@ dipole
use dipole formulation for DAS
@ strainIntegration
use strain integration for DAS