GEOS
CompositionalMultiphaseWellKernels.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_FLUIDFLOW_WELLS_COMPOSITIONALMULTIPHASEWELLKERNELS_HPP
21 #define GEOS_PHYSICSSOLVERS_FLUIDFLOW_WELLS_COMPOSITIONALMULTIPHASEWELLKERNELS_HPP
22 
23 #include "codingUtilities/Utilities.hpp"
24 #include "common/DataTypes.hpp"
25 #include "common/GEOS_RAJA_Interface.hpp"
26 #include "constitutive/fluid/multifluid/MultiFluidBase.hpp"
27 #include "constitutive/fluid/multifluid/MultiFluidFields.hpp"
28 #include "constitutive/relativePermeability/RelativePermeabilityBase.hpp"
29 #include "constitutive/relativePermeability/RelativePermeabilityFields.hpp"
32 #include "mesh/WellElementSubRegion.hpp"
41 #include "physicsSolvers/fluidFlow/wells/WellControls.hpp"
43 
44 namespace geos
45 {
46 
47 
48 namespace compositionalMultiphaseWellKernels
49 {
50 
51 static constexpr real64 minDensForDivision = 1e-10;
52 
53 // tag to access well and reservoir elements in perforation rates computation
55 {
56  static constexpr integer RES = 0;
57  static constexpr integer WELL = 1;
58 };
59 
60 // tag to access the next and current well elements of a connection
61 struct ElemTag
62 {
63  static constexpr integer CURRENT = 0;
64  static constexpr integer NEXT = 1;
65 };
66 
67 // define the column offset of the derivatives
68 struct ColOffset
69 {
70  static constexpr integer DPRES = 0;
71  static constexpr integer DCOMP = 1;
72 };
73 
74 template< integer NC, integer IS_THERMAL >
76 
77 template< integer NC >
78 struct ColOffset_WellJac< NC, 0 >
79 {
80  static constexpr integer dP = 0;
81  static constexpr integer dC = 1;
82  static constexpr integer dQ = dC + NC;
83  static integer constexpr nDer = dQ + 1;
84 
85 };
86 
87 template< integer NC >
88 struct ColOffset_WellJac< NC, 1 >
89 {
90  static constexpr integer dP = 0;
91  static constexpr integer dC = 1;
92  static constexpr integer dQ = dC + NC;
93  static constexpr integer dT = dQ+1;
95  static integer constexpr nDer = dT + 1;
96 };
97 
98 // define the row offset of the residual equations
99 struct RowOffset
100 {
101  static constexpr integer CONTROL = 0;
102  static constexpr integer MASSBAL = 1;
103 };
104 
105 template< integer NC, integer IS_THERMAL >
107 
108 template< integer NC >
109 struct RowOffset_WellJac< NC, 0 >
110 {
111  static constexpr integer CONTROL = 0;
112  static constexpr integer MASSBAL = 1;
113  static constexpr integer VOLBAL = MASSBAL + NC;
114  static constexpr integer nEqn = VOLBAL+1;
115 };
116 
117 template< integer NC >
118 struct RowOffset_WellJac< NC, 1 >
119 {
120  static constexpr integer CONTROL = 0;
121  static constexpr integer MASSBAL = 1;
122  static constexpr integer VOLBAL = MASSBAL + NC;
123  static constexpr integer ENERGYBAL = VOLBAL+1;
124  static constexpr integer nEqn = ENERGYBAL+1;
125 
126 };
127 /******************************** ControlEquationHelper ********************************/
129 {
132 
134  inline
135  static
136  void
137  switchControl( bool const isProducer,
138  WellControls::Control const & inputControl,
139  WellControls::Control const & currentControl,
140  integer const phasePhaseIndex,
141  real64 const & targetBHP,
142  real64 const & targetPhaseRate,
143  real64 const & targetTotalRate,
144  real64 const & targetMassRate,
145  real64 const & currentBHP,
146  arrayView1d< real64 const > const & currentPhaseVolRate,
147  real64 const & currentTotalVolRate,
148  real64 const & currentMassRate,
149  WellControls::Control & newControl );
150 
151  template< integer NC, integer IS_THERMAL >
153  inline
154  static void
155  compute( globalIndex const rankOffset,
156  WellControls::Control const currentControl,
157  integer const targetPhaseIndex,
158  real64 const & targetBHP,
159  real64 const & targetPhaseRate,
160  real64 const & targetTotalRate,
161  real64 const & targetMassRate,
162  real64 const & currentBHP,
163  arrayView1d< real64 const > const & dCurrentBHP,
164  arrayView1d< real64 const > const & currentPhaseVolRate,
165  arrayView2d< real64 const > const & dCurrentPhaseVolRate,
166  real64 const & currentTotalVolRate,
167  arrayView1d< real64 const > const & dCurrentTotalVolRate,
168  real64 const & massDensity,
169  globalIndex const dofNumber,
170  CRSMatrixView< real64, globalIndex const > const & localMatrix,
171  arrayView1d< real64 > const & localRhs );
172 
173 };
174 
175 /******************************** PressureRelationKernel ********************************/
176 
178 {
179  using Deriv = constitutive::multifluid::DerivativeOffset;
183 
184  template< integer NC, integer IS_THERMAL >
186  inline
187  static void
188  compute( real64 const & gravCoef,
189  real64 const & gravCoefNext,
190  real64 const & pres,
191  real64 const & presNext,
192  real64 const & totalMassDens,
193  real64 const & totalMassDensNext,
196  real64 & localPresRel,
197  real64 ( &localPresRelJacobian )[2*(NC+1+IS_THERMAL)] );
198 
199  template< integer NC, integer IS_THERMAL >
200  static void
201  launch( localIndex const size,
202  globalIndex const rankOffset,
203  bool const isLocallyOwned,
204  localIndex const iwelemControl,
205  integer const targetPhaseIndex,
206  WellControls const & wellControls,
207  real64 const & time,
208  arrayView1d< integer const > const elemStatus,
209  arrayView1d< globalIndex const > const & wellElemDofNumber,
210  arrayView1d< real64 const > const & wellElemGravCoef,
211  arrayView1d< localIndex const > const & nextWellElemIndex,
212  arrayView1d< real64 const > const & wellElemPressure,
213  arrayView1d< real64 const > const & wellElemTotalMassDens,
214  arrayView2d< real64 const, compflow::USD_FLUID_DC > const & dWellElemTotalMassDens,
215  bool & controlHasSwitched,
216  CRSMatrixView< real64, globalIndex const > const & localMatrix,
217  arrayView1d< real64 > const & localRhs );
218 
219 };
220 
221 /******************************** VolumeBalanceKernel ********************************/
222 
224 {
225 
228 
229  template< integer NC >
231  inline
232  static void
233  compute( integer const numPhases,
234  real64 const & volume,
237  real64 & localVolBalance,
238  real64 ( &localVolBalanceJacobian )[NC+1] );
239 
240  template< integer NC >
241  static void
242  launch( localIndex const size,
243  integer const numPhases,
244  globalIndex const rankOffset,
245  arrayView1d< globalIndex const > const & wellElemDofNumber,
246  arrayView1d< integer const > const & wellElemGhostRank,
247  arrayView2d< real64 const, compflow::USD_PHASE > const & wellElemPhaseVolFrac,
248  arrayView3d< real64 const, compflow::USD_PHASE_DC > const & dWellElemPhaseVolFrac,
249  arrayView1d< real64 const > const & wellElemVolume,
250  CRSMatrixView< real64, globalIndex const > const & localMatrix,
251  arrayView1d< real64 > const & localRhs );
252 
253 };
254 
255 /******************************** PresTempCompFracInitializationKernel ********************************/
256 
258 {
259 
260  using CompFlowAccessors =
261  StencilAccessors< fields::flow::pressure,
262  fields::flow::temperature,
263  fields::flow::globalCompDensity,
264  fields::flow::phaseVolumeFraction >;
265 
266  using MultiFluidAccessors =
267  StencilMaterialAccessors< constitutive::MultiFluidBase,
268  fields::multifluid::phaseMassDensity >;
269 
270 
278  template< typename VIEWTYPE >
280 
281  static void
282  launch( localIndex const perforationSize,
283  localIndex const subRegionSize,
284  integer const numComponents,
285  integer const numPhases,
286  WellControls const & wellControls,
287  real64 const & currentTime,
293  arrayView1d< localIndex const > const & resElementRegion,
294  arrayView1d< localIndex const > const & resElementSubRegion,
295  arrayView1d< localIndex const > const & resElementIndex,
296  arrayView1d< real64 const > const & perfGravCoef,
297  arrayView1d< integer const > const & perfState,
298  arrayView1d< real64 const > const & wellElemGravCoef,
299  arrayView1d< real64 > const & wellElemPres,
300  arrayView1d< real64 > const & wellElemTemp,
301  arrayView2d< real64, compflow::USD_COMP > const & wellElemCompFrac );
302 
303 };
304 
305 /******************************** CompDensInitializationKernel ********************************/
306 
308 {
309 
310  static void
311  launch( localIndex const subRegionSize,
312  integer const numComponents,
313  arrayView2d< real64 const, compflow::USD_COMP > const & wellElemCompFrac,
315  arrayView2d< real64, compflow::USD_COMP > const & wellElemCompDens );
316 
317 };
318 
319 /******************************** RateInitializationKernel ********************************/
320 
322 {
323 
324  static void
325  launch( localIndex const subRegionSize,
326  integer const targetPhaseIndex,
327  WellControls const & wellControls,
328  real64 const & currentTime,
331  arrayView1d< real64 > const & connRate );
332 
333 };
334 
335 
336 /******************************** TotalMassDensityKernel ****************************/
337 
344 template< integer NUM_COMP, integer NUM_PHASE >
346 {
347 public:
348 
350  using Base::numComp;
351 
353  static constexpr integer numPhase = NUM_PHASE;
354 
361  constitutive::MultiFluidBase const & fluid )
362  : Base(),
363  m_phaseVolFrac( subRegion.getField< fields::well::phaseVolumeFraction >() ),
364  m_dPhaseVolFrac( subRegion.getField< fields::well::dPhaseVolumeFraction >() ),
365  m_dCompFrac_dCompDens( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
366  m_phaseMassDens( fluid.phaseMassDensity() ),
367  m_dPhaseMassDens( fluid.dPhaseMassDensity() ),
368  m_totalMassDens( subRegion.getField< fields::well::totalMassDensity >() ),
369  m_dTotalMassDens( subRegion.getField< fields::well::dTotalMassDensity >() )
370  {}
371 
378  template< typename FUNC = NoOpFunc >
380  inline
381  void compute( localIndex const ei,
382  FUNC && totalMassDensityKernelOp = NoOpFunc{} ) const
383  {
384  using Deriv = constitutive::multifluid::DerivativeOffset;
385 
386  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
387  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
388  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
389  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseMassDens = m_phaseMassDens[ei][0];
390  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseMassDens = m_dPhaseMassDens[ei][0];
391  real64 & totalMassDens = m_totalMassDens[ei];
392  arraySlice1d< real64, compflow::USD_FLUID_DC - 1 > dTotalMassDens = m_dTotalMassDens[ei];
393 
394  real64 dMassDens_dC[numComp]{};
395 
396  totalMassDens = 0.0;
397 
398  dTotalMassDens[Deriv::dP]=0.0;
399  for( integer ic = 0; ic < numComp; ++ic )
400  {
401  dTotalMassDens[Deriv::dC+ic]=0.0;
402  }
403 
404  for( integer ip = 0; ip < numPhase; ++ip )
405  {
406  totalMassDens += phaseVolFrac[ip] * phaseMassDens[ip];
407  dTotalMassDens[Deriv::dP] += dPhaseVolFrac[ip][Deriv::dP] * phaseMassDens[ip] + phaseVolFrac[ip] * dPhaseMassDens[ip][Deriv::dP];
408 
409  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseMassDens[ip], dMassDens_dC, Deriv::dC );
410  for( integer ic = 0; ic < numComp; ++ic )
411  {
412  dTotalMassDens[Deriv::dC+ic] += dPhaseVolFrac[ip][Deriv::dC+ic] * phaseMassDens[ip]
413  + phaseVolFrac[ip] * dMassDens_dC[ic];
414  }
415 
416  totalMassDensityKernelOp( ip ); //, phaseVolFrac, dTotalMassDens_dPres, dTotalMassDens_dCompDens );
417  }
418 
419  }
420 
421 protected:
422 
423  // inputs
424 
429 
433 
434  // outputs
435 
439 
440 
441 };
442 
447 {
448 public:
449 
458  template< typename POLICY >
459  static void
460  createAndLaunch( integer const numComp,
461  integer const numPhase,
462  ObjectManagerBase & subRegion,
463  constitutive::MultiFluidBase const & fluid )
464  {
465  if( numPhase == 2 )
466  {
467  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
468  {
469  integer constexpr NUM_COMP = NC();
470  TotalMassDensityKernel< NUM_COMP, 2 > kernel( subRegion, fluid );
471  TotalMassDensityKernel< NUM_COMP, 2 >::template launch< POLICY >( subRegion.size(), kernel );
472  } );
473  }
474  else if( numPhase == 3 )
475  {
476  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
477  {
478  integer constexpr NUM_COMP = NC();
479  TotalMassDensityKernel< NUM_COMP, 3 > kernel( subRegion, fluid );
480  TotalMassDensityKernel< NUM_COMP, 3 >::template launch< POLICY >( subRegion.size(), kernel );
481  } );
482  }
483  }
484 };
485 
486 
487 /******************************** ResidualNormKernel ********************************/
488 
493 {
494 public:
495 
497  using Base::m_minNormalizer;
498  using Base::m_rankOffset;
499  using Base::m_localResidual;
500  using Base::m_dofNumber;
501 
502  ResidualNormKernel( globalIndex const rankOffset,
503  arrayView1d< real64 const > const & localResidual,
504  arrayView1d< globalIndex const > const & dofNumber,
506  integer const numComp,
507  integer const numDof,
508  integer const targetPhaseIndex,
509  WellElementSubRegion const & subRegion,
510  constitutive::MultiFluidBase const & fluid,
511  WellControls const & wellControls,
512  real64 const time,
513  real64 const dt,
514  real64 const minNormalizer )
515  : Base( rankOffset,
516  localResidual,
517  dofNumber,
518  ghostRank,
519  minNormalizer ),
520  m_numComp( numComp ),
521  m_numDof( numDof ),
522  m_targetPhaseIndex( targetPhaseIndex ),
523  m_dt( dt ),
524  m_isLocallyOwned( subRegion.isLocallyOwned() ),
526  m_isProducer( wellControls.isProducer() ),
527  m_currentControl( wellControls.getControl() ),
528  m_targetBHP( wellControls.getTargetBHP( time ) ),
529  m_targetTotalRate( wellControls.getTargetTotalRate( time ) ),
530  m_targetPhaseRate( wellControls.getTargetPhaseRate( time ) ),
531  m_targetMassRate( wellControls.getTargetMassRate( time ) ),
532  m_volume( subRegion.getElementVolume() ),
533  m_phaseDens_n( fluid.phaseDensity_n() ),
534  m_totalDens_n( fluid.totalDensity_n() )
535  {}
536 
538  virtual void computeLinf( localIndex const iwelem,
539  LinfStackVariables & stack ) const override
540  {
542 
543  real64 normalizer = 0.0;
544  for( integer idof = 0; idof < m_numDof; ++idof )
545  {
546 
547  // Step 1: compute a normalizer for the control or pressure equation
548 
549  // for the control equation, we distinguish two cases
550  if( idof == ROFFSET::CONTROL )
551  {
552 
553  // for the top well element, normalize using the current control
554  if( m_isLocallyOwned && iwelem == m_iwelemControl )
555  {
557  {
558  // the residual entry is in pressure units
559  normalizer = m_targetBHP;
560  }
562  {
563  // the residual entry is in volume / time units
564  normalizer = LvArray::math::max( LvArray::math::abs( m_targetTotalRate ), m_minNormalizer );
565  }
567  {
568  // the residual entry is in volume / time units
569  normalizer = LvArray::math::max( LvArray::math::abs( m_targetPhaseRate ), m_minNormalizer );
570  }
572  {
573  // the residual entry is in volume / time units
574  normalizer = LvArray::math::max( LvArray::math::abs( m_targetMassRate ), m_minNormalizer );
575  }
576  }
577  // for the pressure difference equation, always normalize by the BHP
578  else
579  {
580  normalizer = m_targetBHP;
581  }
582  }
583  // Step 2: compute a normalizer for the mass balance equations
584  else if( idof >= ROFFSET::MASSBAL && idof < ROFFSET::MASSBAL + m_numComp )
585  {
586  if( m_isProducer ) // only PHASEVOLRATE is supported for now
587  {
588  // the residual is in mass units
589  normalizer = m_dt * LvArray::math::abs( m_targetPhaseRate ) * m_phaseDens_n[iwelem][0][m_targetPhaseIndex];
590  }
591  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
592  {
594  {
595  normalizer = m_dt * LvArray::math::abs( m_targetMassRate );
596  }
597  else
598  {
599  // the residual is in mass units
600  normalizer = m_dt * LvArray::math::abs( m_targetTotalRate ) * m_totalDens_n[iwelem][0];
601  }
602 
603  }
604 
605  // to make sure that everything still works well if the rate is zero, we add this check
606  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] * m_totalDens_n[iwelem][0] );
607  }
608  // Step 3: compute a normalizer for the volume balance equations
609  else if( idof == ROFFSET::MASSBAL + m_numComp )
610  {
611  if( m_isProducer ) // only PHASEVOLRATE is supported for now
612  {
613  // the residual is in volume units
614  normalizer = m_dt * LvArray::math::abs( m_targetPhaseRate );
615  }
616  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
617  {
619  {
620  normalizer = m_dt * LvArray::math::abs( m_targetMassRate/ m_totalDens_n[iwelem][0] );
621  }
622  else
623  {
624  normalizer = m_dt * LvArray::math::abs( m_targetTotalRate );
625  }
626 
627  }
628 
629  }
630 
631  // to make sure that everything still works well if the rate is zero, we add this check
632  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] );
633 
634  // Step 4: compute the contribution to the residual
635  real64 const val = LvArray::math::abs( m_localResidual[stack.localRow + idof] ) / normalizer;
636  if( val > stack.localValue[0] )
637  {
638  stack.localValue[0] = val;
639  }
640  }
641  }
642 
644  virtual void computeL2( localIndex const iwelem,
645  L2StackVariables & stack ) const override
646  {
647  GEOS_UNUSED_VAR( iwelem, stack );
648  GEOS_ERROR( "The L2 norm is not implemented for CompositionalMultiphaseWell" );
649  }
650 
651 
652 protected:
653 
656 
659 
662 
664  real64 const m_dt;
665 
667  bool const m_isLocallyOwned;
668 
671 
673  bool const m_isProducer;
674 
677  real64 const m_targetBHP;
678  real64 const m_targetTotalRate;
679  real64 const m_targetPhaseRate;
680  real64 const m_targetMassRate;
681 
684 
688 
689 };
690 
695 {
696 public:
697 
714  template< typename POLICY >
715  static void
716  createAndLaunch( integer const numComp,
717  integer const numDof,
718  integer const targetPhaseIndex,
719  globalIndex const rankOffset,
720  string const & dofKey,
721  arrayView1d< real64 const > const & localResidual,
722  WellElementSubRegion const & subRegion,
723  constitutive::MultiFluidBase const & fluid,
724  WellControls const & wellControls,
725  real64 const time,
726  real64 const dt,
727  real64 const minNormalizer,
728  real64 (& residualNorm)[1] )
729  {
730  arrayView1d< globalIndex const > const dofNumber = subRegion.getReference< array1d< globalIndex > >( dofKey );
731  arrayView1d< integer const > const ghostRank = subRegion.ghostRank();
732 
733  ResidualNormKernel kernel( rankOffset, localResidual, dofNumber, ghostRank,
734  numComp, numDof, targetPhaseIndex, subRegion, fluid, wellControls, time, dt, minNormalizer );
735  ResidualNormKernel::launchLinf< POLICY >( subRegion.size(), kernel, residualNorm );
736  }
737 
738 };
739 
740 /******************************** SolutionScalingKernel ********************************/
741 
746 {
747 public:
748 
762  template< typename POLICY >
764  createAndLaunch( real64 const maxRelativePresChange,
765  real64 const maxAbsolutePresChange,
766  real64 const maxCompFracChange,
767  real64 const maxRelativeCompDensChange,
768  globalIndex const rankOffset,
769  integer const numComp,
770  string const dofKey,
771  ElementSubRegionBase & subRegion,
772  arrayView1d< real64 const > const localSolution )
773  {
774  arrayView1d< real64 const > const pressure =
775  subRegion.getField< fields::well::pressure >();
777  subRegion.getField< fields::well::globalCompDensity >();
778  arrayView1d< real64 > pressureScalingFactor =
779  subRegion.getField< fields::well::pressureScalingFactor >();
780  arrayView1d< real64 > compDensScalingFactor =
781  subRegion.getField< fields::well::globalCompDensityScalingFactor >();
782  isothermalCompositionalMultiphaseBaseKernels::
783  SolutionScalingKernel kernel( maxRelativePresChange, maxAbsolutePresChange, maxCompFracChange, maxRelativeCompDensChange, rankOffset,
784  numComp, dofKey, subRegion, localSolution, pressure, compDens, pressureScalingFactor, compDensScalingFactor );
785  return isothermalCompositionalMultiphaseBaseKernels::
786  SolutionScalingKernel::
787  launch< POLICY >( subRegion.size(), kernel );
788  }
789 
790 };
791 
792 /******************************** ElementBasedAssemblyKernel ********************************/
793 
800 template< integer NUM_COMP, integer IS_THERMAL >
802 {
803 public:
806 
807  // Well jacobian column and row indicies
808  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NUM_COMP, IS_THERMAL >;
812  static constexpr integer numComp = NUM_COMP;
813 
815  static constexpr integer numDof = NUM_COMP + 1 + IS_THERMAL;
816 
818  static constexpr integer numEqn = NUM_COMP + 1 + IS_THERMAL;
819 
820 
833  integer const isProducer,
834  globalIndex const rankOffset,
835  string const dofKey,
836  WellElementSubRegion const & subRegion,
837  constitutive::MultiFluidBase const & fluid,
838  CRSMatrixView< real64, globalIndex const > const & localMatrix,
839  arrayView1d< real64 > const & localRhs,
840  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const kernelFlags )
841  : m_numPhases( numPhases ),
842  m_isProducer( isProducer ),
843  m_rankOffset( rankOffset ),
844  m_iwelemControl( subRegion.getTopWellElementIndex() ),
845  m_dofNumber( subRegion.getReference< array1d< globalIndex > >( dofKey ) ),
846  m_elemGhostRank( subRegion.ghostRank() ),
847  m_elemStatus( subRegion.getLocalWellElementStatus() ),
848  m_volume( subRegion.getElementVolume() ),
849  m_dCompFrac_dCompDens( subRegion.getField< fields::flow::dGlobalCompFraction_dGlobalCompDensity >() ),
850  m_phaseVolFrac_n( subRegion.getField< fields::flow::phaseVolumeFraction_n >() ),
851  m_phaseVolFrac( subRegion.getField< fields::flow::phaseVolumeFraction >() ),
852  m_dPhaseVolFrac( subRegion.getField< fields::flow::dPhaseVolumeFraction >() ),
853  m_phaseDens_n( fluid.phaseDensity_n() ),
854  m_phaseDens( fluid.phaseDensity() ),
855  m_dPhaseDens( fluid.dPhaseDensity() ),
856  m_phaseCompFrac_n( fluid.phaseCompFraction_n() ),
857  m_phaseCompFrac( fluid.phaseCompFraction() ),
858  m_dPhaseCompFrac( fluid.dPhaseCompFraction() ),
859  m_compDens( subRegion.getField< fields::flow::globalCompDensity >() ),
860  m_compDens_n( subRegion.getField< fields::flow::globalCompDensity_n >() ),
861  m_localMatrix( localMatrix ),
862  m_localRhs( localRhs ),
863  m_kernelFlags( kernelFlags )
864  {}
865 
871  {
872 public:
873 
874  // volume information (used by both accumulation and volume balance)
875  real64 volume = 0.0;
876 
877  // Residual information
878 
881 
883  globalIndex dofIndices[numDof]{}; // NC compdens + P + thermal
884  globalIndex eqnRowIndices[numDof]{};
885  globalIndex dofColIndices[numDof]{};
886 
889 
892 
893  };
894 
901  bool skipElement( localIndex const ei ) const
902  {
903  return ( m_elemGhostRank( ei ) >= 0 || m_elemStatus[ei]==WellElementSubRegion::WellElemStatus::CLOSED );
904  }
905 
912  void setup( localIndex const ei,
913  StackVariables & stack ) const
914  {
915  // initialize the volume
916  stack.volume = m_volume[ei];
917 
918  // Note row/col indices needed to be consistent with layout of stack.localJacobian
919  // Setup row equation indices for this element ( mass + vol + thermal if valid)
920 
921  // 1) Mass Balance
922  for( integer ic = 0; ic < numComp; ++ic )
923  {
924  stack.eqnRowIndices[ic] = m_dofNumber[ei] + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
925  }
926 // 2) Volume Balance
927  stack.eqnRowIndices[numComp] = m_dofNumber[ei] + WJ_ROFFSET::VOLBAL - m_rankOffset;
928  // 3) Energy Balance
929  if constexpr ( IS_THERMAL )
930  {
931  stack.eqnRowIndices[numComp+1] = m_dofNumber[ei] + WJ_ROFFSET::ENERGYBAL - m_rankOffset;
932  }
933  // Setup equation column indices for this element ( P + COMPDENS + THERMAL if valid)
934  stack.dofColIndices[0] = m_dofNumber[ei] + WJ_COFFSET::dP;
935  for( integer ic = 0; ic < numComp; ++ic )
936  {
937  stack.dofColIndices[ic+1] = m_dofNumber[ei] + WJ_COFFSET::dC+ic;
938  }
939  if constexpr ( IS_THERMAL )
940  {
941  stack.dofColIndices[numComp+1] = m_dofNumber[ei] + WJ_COFFSET::dT;
942  }
943  for( integer jc = 0; jc < numEqn; ++jc )
944  {
945  stack.localResidual[jc] = 0.0;
946  for( integer ic = 0; ic < numDof; ++ic )
947  {
948  stack.localJacobian[jc][ic] = 0.0;
949  }
950  }
951  }
959  template< typename FUNC = NoOpFunc >
962  StackVariables & stack,
963  FUNC && phaseAmountKernelOp = NoOpFunc{} ) const
964  {
965 
966  using Deriv = constitutive::multifluid::DerivativeOffset;
967 
968  // construct the slices for variables accessed multiple times
969  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
970 
971  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac_n = m_phaseVolFrac_n[ei];
972  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
973  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
974 
975  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens_n = m_phaseDens_n[ei][0];
976  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens = m_phaseDens[ei][0];
977  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseDens = m_dPhaseDens[ei][0];
978 
979  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac_n = m_phaseCompFrac_n[ei][0];
980  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac = m_phaseCompFrac[ei][0];
981  arraySlice3d< real64 const, constitutive::multifluid::USD_PHASE_COMP_DC - 2 > dPhaseCompFrac = m_dPhaseCompFrac[ei][0];
982 
983  // temporary work arrays
984  real64 dPhaseAmount[FLUID_PROP_COFFSET::nDer]{};
985  real64 dPhaseAmount_dC[numComp]{};
986  real64 dPhaseCompFrac_dC[numComp]{};
987 
988  // sum contributions to component accumulation from each phase
989  for( integer ip = 0; ip < m_numPhases; ++ip )
990  {
991  real64 const phaseAmount = stack.volume * phaseVolFrac[ip] * phaseDens[ip];
992  real64 const phaseAmount_n = stack.volume * phaseVolFrac_n[ip] * phaseDens_n[ip];
993 
994  dPhaseAmount[FLUID_PROP_COFFSET::dP]=stack.volume * ( dPhaseVolFrac[ip][Deriv::dP] * phaseDens[ip]
995  + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dP] );
996 
997  // assemble density dependence
998  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], dPhaseAmount_dC, Deriv::dC );
999  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], &dPhaseAmount[FLUID_PROP_COFFSET::dC], Deriv::dC );
1000  for( integer jc = 0; jc < numComp; ++jc )
1001  {
1002  dPhaseAmount_dC[jc] = dPhaseAmount_dC[jc] * phaseVolFrac[ip]
1003  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1004  dPhaseAmount_dC[jc] *= stack.volume;
1005  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] = dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] * phaseVolFrac[ip]
1006  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1007  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] *= stack.volume;
1008  }
1009  // ic - index of component whose conservation equation is assembled
1010  // (i.e. row number in local matrix)
1011  for( integer ic = 0; ic < numComp; ++ic )
1012  {
1013  real64 const phaseCompAmount = phaseAmount * phaseCompFrac[ip][ic];
1014  real64 const phaseCompAmount_n = phaseAmount_n * phaseCompFrac_n[ip][ic];
1015 
1016  real64 const dPhaseCompAmount_dP = dPhaseAmount[FLUID_PROP_COFFSET::dP] * phaseCompFrac[ip][ic]
1017  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dP];
1018 
1019  stack.localResidual[ic] += phaseCompAmount - phaseCompAmount_n;
1020  stack.localJacobian[ic][0] += dPhaseCompAmount_dP;
1021 
1022  // jc - index of component w.r.t. whose compositional var the derivative is being taken
1023  // (i.e. col number in local matrix)
1024 
1025  // assemble phase composition dependence
1026  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseCompFrac[ip][ic], dPhaseCompFrac_dC, Deriv::dC );
1027  for( integer jc = 0; jc < numComp; ++jc )
1028  {
1029  real64 const dPhaseCompAmount_dC = dPhaseCompFrac_dC[jc] * phaseAmount
1030  + phaseCompFrac[ip][ic] * dPhaseAmount[FLUID_PROP_COFFSET::dC+jc];
1031 
1032  stack.localJacobian[ic][jc + 1] += dPhaseCompAmount_dC;
1033  }
1034  }
1035  if constexpr ( IS_THERMAL )
1036  {
1037  dPhaseAmount[FLUID_PROP_COFFSET::dT] = stack.volume * (dPhaseVolFrac[ip][Deriv::dT] * phaseDens[ip] + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dT] );
1038  for( integer ic = 0; ic < numComp; ++ic )
1039  {
1040  // assemble the derivatives of the component mass balance equations with respect to temperature
1041  stack.localJacobian[ic][numComp+1] += dPhaseAmount[FLUID_PROP_COFFSET::dT] * phaseCompFrac[ip][ic]
1042  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dT];
1043  }
1044  }
1045  // call the lambda in the phase loop to allow the reuse of the phase amounts and their derivatives
1046  // possible use: assemble accumulation term of the energy equation for this phase
1047  phaseAmountKernelOp( ip, phaseAmount, phaseAmount_n, dPhaseAmount );
1048 
1049  }
1050 
1051  // check zero diagonal (works only in debug)
1052  /*
1053  for( integer ic = 0; ic < numComp; ++ic )
1054  {
1055  GEOS_ASSERT_MSG ( LvArray::math::abs( stack.localJacobian[ic][ic] ) > minDensForDivision,
1056  GEOS_FMT( "Zero diagonal in Jacobian: equation {}, value = {}", ic, stack.localJacobian[ic][ic] ) );
1057  }
1058  */
1059  }
1060 
1061 
1072  StackVariables & stack ) const
1073  {
1074  using Deriv = constitutive::multifluid::DerivativeOffset;
1075 
1076  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
1077  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
1078 
1079  real64 oneMinusPhaseVolFracSum = 1.0;
1080 
1081  // sum contributions to component accumulation from each phase
1082 // Note localJacobian stores equation balances in order of component/vol/enerqy
1083  // These are mapped to solver orderings with indicies setup in stack variables
1084  for( integer ip = 0; ip < m_numPhases; ++ip )
1085  {
1086  oneMinusPhaseVolFracSum -= phaseVolFrac[ip];
1087  stack.localJacobian[numComp][0] -= dPhaseVolFrac[ip][Deriv::dP];
1088 
1089  for( integer jc = 0; jc < numComp; ++jc )
1090  {
1091  stack.localJacobian[numComp][jc+1] -= dPhaseVolFrac[ip][Deriv::dC+jc];
1092  }
1093 
1094  if constexpr ( IS_THERMAL)
1095  {
1096  stack.localJacobian[numComp][numComp+1] -= dPhaseVolFrac[ip][Deriv::dT];
1097  }
1098 
1099  }
1100  // scale saturation-based volume balance by pore volume (for better scaling w.r.t. other equations)
1101  stack.localResidual[numComp] = stack.volume * oneMinusPhaseVolFracSum;
1102  for( integer idof = 0; idof < numComp+1+IS_THERMAL; ++idof )
1103  {
1104  stack.localJacobian[numComp][idof] *= stack.volume;
1105  }
1106  }
1107 
1114  void complete( localIndex const ei, //GEOS_UNUSED_PARAM( ei ),
1115  StackVariables & stack ) const
1116  {
1117  using namespace compositionalMultiphaseUtilities;
1118 
1119  integer const numRows = numComp+1+ IS_THERMAL;
1120 
1121  if constexpr ( IS_THERMAL)
1122  {
1123  if( ei == m_iwelemControl && !m_isProducer )
1124  {
1125  // For top segment energy balance eqn replaced with T(n+1) - T = 0
1126  // No other energy balance derivatives
1127  // Assumption is global index == 0 is top segment with fixed temp BC
1128 
1129  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1130  {
1131  stack.localJacobian[numRows-1][i] = 0.0;
1132  }
1133  // constant Temperature
1134  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1135  stack.localJacobian[i][numRows-1] = 0.0;
1136  stack.localJacobian[numRows-1][numRows-1] = 1.0;
1137 
1138  stack.localResidual[numRows-1]=0.0;
1139  }
1140  }
1141 
1142  if( m_kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) )
1143  {
1144  // apply equation/variable change transformation to the component mass balance equations
1145  real64 work[numComp + 1 + IS_THERMAL]{};
1146  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numComp, numComp+1+ IS_THERMAL, stack.localJacobian, work );
1147  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numComp, stack.localResidual );
1148  }
1149 
1150  // add contribution to residual and jacobian into:
1151  // - the component mass balance equations (i = 0 to i = numComp-1)
1152  // - the volume balance equations (i = numComp)
1153  // note that numDof includes derivatives wrt temperature if this class is derived in ThermalKernels
1154 
1155  for( integer i = 0; i < numRows; ++i )
1156  {
1157  m_localRhs[stack.eqnRowIndices[i]] += stack.localResidual[i];
1158  m_localMatrix.addToRow< serialAtomic >( stack.eqnRowIndices[i],
1159  stack.dofColIndices,
1160  stack.localJacobian[i],
1161  numComp+1+ IS_THERMAL );
1162  }
1163 
1164  }
1165 
1173  template< typename POLICY, typename KERNEL_TYPE >
1174  static void
1175  launch( localIndex const numElems,
1176  KERNEL_TYPE const & kernelComponent )
1177  {
1179  forAll< POLICY >( numElems, [=] GEOS_HOST_DEVICE ( localIndex const ei )
1180  {
1181  if( kernelComponent.skipElement( ei ) )
1182  {
1183  return;
1184  }
1185 
1186  typename KERNEL_TYPE::StackVariables stack;
1187 
1188  kernelComponent.setup( ei, stack );
1189  kernelComponent.computeAccumulation( ei, stack );
1190  kernelComponent.computeVolumeBalance( ei, stack );
1191  kernelComponent.complete( ei, stack );
1192  } );
1193  }
1194 
1195 protected:
1196 
1199 
1202 
1205 
1208 
1211 
1214 
1217 
1220 
1223  arrayView2d< real64 const > const m_porosity;
1224  arrayView2d< real64 const > const m_dPoro_dPres;
1225 
1228 
1233 
1238 
1243 
1244  // Views on component densities
1247 
1252 
1253  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const m_kernelFlags;
1254 };
1255 
1256 
1261 {
1262 public:
1275  template< typename POLICY >
1276  static void
1277  createAndLaunch( localIndex const numComps,
1278  localIndex const numPhases,
1279  integer const isProducer,
1280  globalIndex const rankOffset,
1281  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags,
1282  string const dofKey,
1283  WellElementSubRegion const & subRegion,
1284  constitutive::MultiFluidBase const & fluid,
1285  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1286  arrayView1d< real64 > const & localRhs )
1287  {
1288  geos::internal::kernelLaunchSelectorCompThermSwitch( numComps, 0, [&]( auto NC, auto IS_THERMAL )
1289  {
1290  localIndex constexpr NUM_COMP = NC();
1291 
1292  integer constexpr istherm = IS_THERMAL();
1293 
1295  kernel( numPhases, isProducer, rankOffset, dofKey, subRegion, fluid, localMatrix, localRhs, kernelFlags );
1297  launch< POLICY, ElementBasedAssemblyKernel< NUM_COMP, istherm > >( subRegion.size(), kernel );
1298  } );
1299  }
1300 };
1309 template< integer NC, integer IS_THERMAL >
1311 {
1312 public:
1313 
1317 
1318  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1321 
1322  using CP_Deriv = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1324  static constexpr integer numComp = NC;
1325 
1327  static constexpr integer numDof = WJ_COFFSET::nDer;
1328 
1330  static constexpr integer numEqn = NC;
1331 
1332  static constexpr integer maxNumElems = 2;
1333  static constexpr integer maxStencilSize = 2;
1349  globalIndex const rankOffset,
1350  string const wellDofKey,
1351  WellControls const & wellControls,
1352  WellElementSubRegion const & subRegion,
1353  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1354  arrayView1d< real64 > const & localRhs,
1355  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags )
1356  :
1357  m_dt( dt ),
1358  m_rankOffset( rankOffset ),
1359  m_wellElemDofNumber ( subRegion.getReference< array1d< globalIndex > >( wellDofKey ) ),
1360  m_nextWellElemIndex ( subRegion.getReference< array1d< localIndex > >( WellElementSubRegion::viewKeyStruct::nextWellElementIndexString()) ),
1361  m_elemStatus( subRegion.getLocalWellElementStatus() ),
1362  m_connRate ( subRegion.getField< fields::well::mixtureConnectionRate >() ),
1363  m_wellElemCompFrac ( subRegion.getField< fields::well::globalCompFraction >() ),
1364  m_dWellElemCompFrac_dCompDens ( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
1365  m_localMatrix( localMatrix ),
1366  m_localRhs ( localRhs ),
1367  m_useTotalMassEquation ( kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) ),
1368  m_isProducer ( wellControls.isProducer() ),
1369  m_injection ( wellControls.getInjectionStream() )
1370  {}
1371 
1373  {
1374 public:
1375 
1383  : stencilSize( size ),
1384  numConnectedElems( 2 ),
1385  dofColIndices( size * numDof )
1386  {}
1387 
1388  // Stencil information
1389  localIndex const stencilSize;
1392 
1393 
1394  // edge indexes
1395  localIndex iwelemUp;
1396  localIndex iwelemNext;
1397  localIndex iwelemCurrent;
1398  globalIndex offsetUp;
1399  globalIndex offsetCurrent;
1400  globalIndex offsetNext;
1401  // Local degrees of freedom and local residual/jacobian
1402 
1405 
1412  };
1413 
1414 
1416  bool skipElement( localIndex const ei ) const
1417  {
1418  return ( m_elemStatus[ei]==WellElementSubRegion::WellElemStatus::CLOSED );
1419  }
1420 
1427  inline
1428  void setup( localIndex const iconn,
1429  StackVariables & stack ) const
1430  {
1431  stack.numConnectedElems=2;
1432  if( m_nextWellElemIndex[iconn] <0 )
1433  {
1434  stack.numConnectedElems = 1;
1435  }
1436 
1437  stack.localFlux.resize( stack.numConnectedElems*numEqn );
1438  stack.localFluxJacobian.resize( stack.numConnectedElems * numEqn, stack.stencilSize * numDof );
1439  stack.localFluxJacobian_dQ.resize( stack.numConnectedElems * numEqn, 1 );
1440 
1441  }
1442 
1449  inline
1450  void complete( localIndex const iconn,
1451  StackVariables & stack ) const
1452  {
1453  GEOS_UNUSED_VAR( iconn );
1454  using namespace compositionalMultiphaseUtilities;
1455  if( stack.numConnectedElems ==1 )
1456  {
1457  // Setup Jacobian global row indicies
1458  // equations for COMPONENT + ENERGY balances
1459  globalIndex oneSidedEqnRowIndices[numEqn]{};
1460  for( integer ic = 0; ic < NC; ++ic )
1461  {
1462  oneSidedEqnRowIndices[ic] = stack.offsetUp + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1463  }
1464 
1465  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1466  globalIndex oneSidedDofColIndices_dPresCompTempUp[CP_Deriv::nDer]{};
1467  globalIndex oneSidedDofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1468  // Note localFluxJacobian cols are stored using CP_Deriv order (dP dC or dP dT dC)
1469  int ioff=0;
1470  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1471 
1472  if constexpr ( IS_THERMAL )
1473  {
1474  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1475  }
1476  for( integer jdof = 0; jdof < NC; ++jdof )
1477  {
1478  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1479  }
1480  if( m_useTotalMassEquation > 0 )
1481  {
1482  // Apply equation/variable change transformation(s)
1483  real64 work[CP_Deriv::nDer]{};
1484  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, 1, stack.localFluxJacobian_dQ, work );
1485  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, CP_Deriv::nDer, stack.localFluxJacobian, work );
1486  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, stack.localFlux );
1487  }
1488  for( integer i = 0; i < numEqn; ++i )
1489  {
1490  if( oneSidedEqnRowIndices[i] >= 0 && oneSidedEqnRowIndices[i] < m_localMatrix.numRows() )
1491  {
1492  m_localMatrix.addToRow< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1493  &oneSidedDofColIndices_dRate,
1494  stack.localFluxJacobian_dQ[i],
1495  1 );
1496  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1497  oneSidedDofColIndices_dPresCompTempUp,
1498  stack.localFluxJacobian[i],
1499  CP_Deriv::nDer );
1500  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[oneSidedEqnRowIndices[i]], stack.localFlux[i] );
1501  }
1502  }
1503  }
1504  else
1505  {
1506  // Setup Jacobian global row indicies
1507  // equations for COMPONENT + ENERGY balances
1508  globalIndex eqnRowIndices[2*numEqn]{};
1509 
1510  for( integer ic = 0; ic < NC; ++ic )
1511  {
1512  // mass balance equations for all components
1513  eqnRowIndices[TAG::NEXT *numEqn+ic] = stack.offsetNext + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1514  eqnRowIndices[TAG::CURRENT *numEqn+ic] = stack.offsetCurrent + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1515  }
1516 
1517  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1518  globalIndex dofColIndices_dPresCompUp[CP_Deriv::nDer]{};
1519  globalIndex dofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1520 
1521  int ioff=0;
1522  // Indice storage order reflects local jac col storage order CP::Deriv order P T DENS
1523  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1524 
1525  if constexpr ( IS_THERMAL )
1526  {
1527  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1528  }
1529  for( integer jdof = 0; jdof < NC; ++jdof )
1530  {
1531  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1532  }
1533 
1534 
1535  if( m_useTotalMassEquation > 0 )
1536  {
1537  // Apply equation/variable change transformation(s)
1538  real64 work[CP_Deriv::nDer]{};
1539  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, 1, 2, stack.localFluxJacobian_dQ, work );
1540  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, CP_Deriv::nDer, 2, stack.localFluxJacobian, work );
1541  shiftBlockElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, numEqn, 2, stack.localFlux );
1542  }
1543  // Note this updates diag and offdiag
1544  for( integer i = 0; i < 2*NC; ++i )
1545  {
1546  if( eqnRowIndices[i] >= 0 && eqnRowIndices[i] < m_localMatrix.numRows() )
1547  {
1548  m_localMatrix.addToRow< parallelDeviceAtomic >( eqnRowIndices[i],
1549  &dofColIndices_dRate,
1550  stack.localFluxJacobian_dQ[i],
1551  1 );
1552  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( eqnRowIndices[i],
1553  dofColIndices_dPresCompUp,
1554  stack.localFluxJacobian[i],
1555  CP_Deriv::nDer );
1556  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[eqnRowIndices[i]], stack.localFlux[i] );
1557 
1558  }
1559  }
1560  }
1561  }
1562 
1564  inline
1565  void
1566  computeExit( real64 const & dt,
1567  real64 const ( &compFlux )[NC ],
1568  StackVariables & stack,
1569  real64 ( & dCompFlux)[NC][numDof] ) const
1570  {
1571  for( integer ic = 0; ic < NC; ++ic )
1572  {
1573  stack.localFlux[ic] = -dt * compFlux[ic];
1574  // derivative with respect to rate
1575  stack.localFluxJacobian_dQ[ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1576  // derivative with respect to upstream pressure
1577  stack.localFluxJacobian[ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1578  // derivatives with respect to upstream component densities
1579  for( integer jdof = 0; jdof < NC; ++jdof )
1580  {
1581  stack.localFluxJacobian[ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1582  }
1583  if constexpr ( IS_THERMAL )
1584  {
1585  stack.localFluxJacobian[ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1586  }
1587  }
1588  }
1589 
1591  inline
1592  void
1593  compute( real64 const & dt,
1594  real64 const ( &compFlux )[NC ],
1595  StackVariables & stack,
1596  real64 ( & dCompFlux)[(NC )][numDof] ) const
1597  {
1598  // flux terms
1599  for( integer ic = 0; ic < NC; ++ic )
1600  {
1601  stack.localFlux[TAG::NEXT * NC +ic] = dt * compFlux[ic];
1602  stack.localFlux[TAG::CURRENT * NC +ic] = -dt * compFlux[ic];
1603  // derivative with respect to rate
1604  stack.localFluxJacobian_dQ[TAG::NEXT * NC+ ic][0] = dt * dCompFlux[ic][WJ_COFFSET::dQ];
1605  stack.localFluxJacobian_dQ[TAG::CURRENT * NC + +ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1606 
1607  // derivative with respect to upstream pressure
1608  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dP] = dt * dCompFlux[ic][WJ_COFFSET::dP];
1609  stack.localFluxJacobian[TAG::CURRENT * NC+ ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1610 
1611  if constexpr ( IS_THERMAL )
1612  {
1613  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dT] = dt * dCompFlux[ic][WJ_COFFSET::dT];
1614  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1615  }
1616 
1617  // derivatives with respect to upstream component densities
1618  for( integer jdof = 0; jdof < NC; ++jdof )
1619  {
1620  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dC+jdof] = dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1621  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1622  }
1623  }
1624  }
1625 
1633  template< typename FUNC = NoOpFunc >
1635  inline
1636  void computeFlux( localIndex const iwelem,
1637  StackVariables & stack,
1638  FUNC && compFluxKernelOp = NoOpFunc{} ) const
1639  {
1640 
1641  using namespace compositionalMultiphaseUtilities;
1642 
1643  // create local work arrays
1644  real64 compFracUp[NC]{};
1645  real64 compFlux[NC]{};
1646  real64 dComp[NC][NC];
1647  real64 dCompFlux[NC][numDof]{};
1648  for( integer ic = 0; ic < NC; ++ic )
1649  {
1650  for( integer jc = 0; jc < NC; ++jc )
1651  {
1652  dComp[ic][jc]=0.0;
1653  }
1654  }
1655  // Step 1) decide the upwind well element
1656 
1657  /* currentConnRate < 0 flow from iwelem to iwelemNext
1658  * currentConnRate > 0 flow from iwelemNext to iwelem
1659  * With this convention, currentConnRate < 0 at the last connection for a producer
1660  * currentConnRate > 0 at the last connection for a injector
1661  */
1662 
1663  localIndex iwelemNext = m_nextWellElemIndex[iwelem];
1664 
1665  // Current element is open and next element is closed => no dependency on next segment
1666  if( iwelemNext >= 0 )
1667  {
1668  if( m_elemStatus[iwelemNext]==WellElementSubRegion::WellElemStatus::CLOSED )
1669  {
1670  iwelemNext = -1;
1671  }
1672  }
1673 
1674  real64 const currentConnRate = m_connRate[iwelem];
1675  localIndex iwelemUp = -1;
1676  if( iwelemNext < 0 && !m_isProducer ) // exit connection, injector
1677  {
1678  // we still need to define iwelemUp for Jacobian assembly
1679  iwelemUp = iwelem;
1680 
1681  // just copy the injection stream into compFrac
1682  for( integer ic = 0; ic < NC; ++ic )
1683  {
1684  compFracUp[ic] = m_injection[ic];
1685  for( integer jc = 0; jc < NC; ++jc )
1686  {
1687  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1688  }
1689  for( integer jc = 0; jc < NC; ++jc )
1690  {
1691  dCompFlux[ic][WJ_COFFSET::dC+jc] = 0.0;
1692  }
1693  }
1694  }
1695  else
1696  {
1697  // first set iwelemUp to the upstream cell
1698  if( ( iwelemNext < 0 && m_isProducer ) // exit connection, producer
1699  || currentConnRate < 0 ) // not an exit connection, iwelem is upstream
1700  {
1701  iwelemUp = iwelem;
1702  }
1703  else // not an exit connection, iwelemNext is upstream
1704  {
1705  iwelemUp = iwelemNext;
1706  }
1707 
1708  // copy the vars of iwelemUp into compFrac
1709  for( integer ic = 0; ic < NC; ++ic )
1710  {
1711  compFracUp[ic] = m_wellElemCompFrac[iwelemUp][ic];
1712  for( integer jc = 0; jc < NC; ++jc )
1713  {
1714  dCompFlux[ic][WJ_COFFSET::dC+jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1715  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1716  }
1717  }
1718  }
1719 
1720  // Step 2) compute upstream transport coefficient
1721 
1722  for( integer ic = 0; ic < NC; ++ic )
1723  {
1724  compFlux[ic] = compFracUp[ic] * currentConnRate;
1725  dCompFlux[ic][WJ_COFFSET::dQ] = compFracUp[ic];
1726  // none of these quantities depend on pressure
1727  dCompFlux[ic][WJ_COFFSET::dP] = 0.0;
1728  if constexpr ( IS_THERMAL )
1729  {
1730  dCompFlux[ic][WJ_COFFSET::dT] = 0.0;
1731  }
1732  for( integer jc = 0; jc < NC; ++jc )
1733  {
1734  dCompFlux[ic][WJ_COFFSET::dC+jc] = dCompFlux[ic][WJ_COFFSET::dC+jc] * currentConnRate;
1735  }
1736  }
1737 
1738  stack.offsetUp = m_wellElemDofNumber[iwelemUp];
1739  stack.iwelemUp = iwelemUp;
1740  stack.offsetCurrent = m_wellElemDofNumber[iwelem];
1741  stack.iwelemCurrent= iwelem;
1742 
1743 
1744  if( iwelemNext < 0 ) // exit connection
1745  {
1746  // for this case, we only need NC mass conservation equations
1747  computeExit ( m_dt,
1748  compFlux,
1749  stack,
1750  dCompFlux );
1751 
1752  }
1753  else // not an exit connection
1754  {
1755  compute( m_dt,
1756  compFlux,
1757  stack,
1758  dCompFlux
1759  );
1760  stack.offsetNext = m_wellElemDofNumber[iwelemNext];
1761  }
1762  stack.iwelemNext = iwelemNext;
1763  compFluxKernelOp( iwelemNext, iwelemUp, currentConnRate, dComp );
1764  }
1765 
1766 
1774  template< typename POLICY, typename KERNEL_TYPE >
1775  static void
1776  launch( localIndex const numElements,
1777  KERNEL_TYPE const & kernelComponent )
1778  {
1780  forAll< POLICY >( numElements, [=] GEOS_HOST_DEVICE ( localIndex const ie )
1781  {
1782  typename KERNEL_TYPE::StackVariables stack( 1 );
1783  if( kernelComponent.skipElement( ie ))
1784  {
1785  return;
1786  }
1787  kernelComponent.setup( ie, stack );
1788  kernelComponent.computeFlux( ie, stack );
1789  kernelComponent.complete( ie, stack );
1790  } );
1791  }
1792 
1793 protected:
1795  real64 const m_dt;
1798 
1803 
1806 
1809 
1810 
1815 
1820 
1823 
1825  bool const m_isProducer;
1826 
1829 
1830 
1831 };
1832 
1837 {
1838 public:
1839 
1853  template< typename POLICY >
1854  static void
1855  createAndLaunch( integer const numComps,
1856  real64 const dt,
1857  globalIndex const rankOffset,
1858  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags,
1859  string const dofKey,
1860  WellControls const & wellControls,
1861  WellElementSubRegion const & subRegion,
1862  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1863  arrayView1d< real64 > const & localRhs )
1864  {
1865  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComps, [&]( auto NC )
1866  {
1867  integer constexpr NUM_COMP = NC();
1868 
1869  using kernelType = FaceBasedAssemblyKernel< NUM_COMP, 0 >;
1870  kernelType kernel( dt, rankOffset, dofKey, wellControls, subRegion, localMatrix, localRhs, kernelFlags );
1871  kernelType::template launch< POLICY >( subRegion.size(), kernel );
1872  } );
1873  }
1874 };
1875 } // end namespace compositionalMultiphaseWellKernels
1876 
1877 } // end namespace geos
1878 
1879 #endif //GEOS_PHYSICSSOLVERS_FLUIDFLOW_WELLS_COMPOSITIONALMULTIPHASEWELLKERNELS_HPP
#define GEOS_HOST_DEVICE
Marks a host-device function.
Definition: GeosxMacros.hpp:49
#define GEOS_UNUSED_VAR(...)
Mark an unused variable and silence compiler warnings.
Definition: GeosxMacros.hpp:84
#define GEOS_ERROR(msg)
Raise a hard error and terminate the program.
Definition: Logger.hpp:157
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
typename ElementViewAccessor< VIEWTYPE >::NestedViewTypeConst ElementViewConst
The ElementViewAccessor at the ElementRegionManager level is the type resulting from ElementViewAcces...
arrayView1d< real64 const > getElementVolume() const
Get the volume of each element in this subregion.
The ObjectManagerBase is the base object of all object managers in the mesh data hierachy.
array1d< integer > const & ghostRank()
Get the ghost information of each object.
GEOS_DECLTYPE_AUTO_RETURN getField() const
Get a view to the field associated with a trait from this ObjectManagerBase.
A struct to automatically construct and store element view accessors.
A struct to automatically construct and store element view accessors.
This class describes the controls used to operate a well.
bool isProducer() const
Is the well a producer?
Control getControl() const
Get the control type for the well.
real64 getTargetPhaseRate(real64 const &currentTime) const
Get the target phase rate.
real64 getTargetTotalRate(real64 const &currentTime) const
Get the target total rate.
real64 getTargetMassRate(real64 const &currentTime) const
Get the target mass rate.
real64 getTargetBHP(real64 const &currentTime) const
Get the target bottom hole pressure value.
This class describes a collection of local well elements and perforations.
bool isLocallyOwned() const
Check if well is owned by current rank.
localIndex getTopWellElementIndex() const
Get for the top element index.
static void createAndLaunch(localIndex const numComps, localIndex const numPhases, integer const isProducer, globalIndex const rankOffset, BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags, string const dofKey, WellElementSubRegion const &subRegion, constitutive::MultiFluidBase const &fluid, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
Create a new kernel and launch.
Define the interface for the assembly kernel in charge of accumulation and volume balance.
static constexpr integer numEqn
Compute time value for the number of equations mass bal + vol bal + energy bal.
ElementBasedAssemblyKernel(localIndex const numPhases, integer const isProducer, globalIndex const rankOffset, string const dofKey, WellElementSubRegion const &subRegion, constitutive::MultiFluidBase const &fluid, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs, BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const kernelFlags)
Constructor.
arrayView1d< globalIndex const > const m_dofNumber
View on the dof numbers.
arrayView3d< real64 const, constitutive::multifluid::USD_PHASE > const m_phaseDens_n
Views on the phase densities.
GEOS_HOST_DEVICE void computeAccumulation(localIndex const ei, StackVariables &stack, FUNC &&phaseAmountKernelOp=NoOpFunc{}) const
Compute the local accumulation contributions to the residual and Jacobian.
GEOS_HOST_DEVICE bool skipElement(localIndex const ei) const
Getter for the element.
arrayView4d< real64 const, constitutive::multifluid::USD_PHASE_COMP > const m_phaseCompFrac_n
Views on the phase component fraction.
arrayView3d< real64 const, compflow::USD_COMP_DC > const m_dCompFrac_dCompDens
Views on the derivatives of comp fractions wrt component density.
arrayView1d< real64 const > const m_volume
View on the element volumes.
GEOS_HOST_DEVICE void complete(localIndex const ei, StackVariables &stack) const
Performs the complete phase for the kernel.
static constexpr integer numComp
Compile time value for the number of components.
static void launch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Performs the kernel launch.
arrayView2d< real64 const, compflow::USD_PHASE > const m_phaseVolFrac_n
Views on the phase volume fractions.
static constexpr integer numDof
Number of Dof's set in this kernal - no dQ in accum.
GEOS_HOST_DEVICE void computeVolumeBalance(localIndex const ei, StackVariables &stack) const
Compute the local volume balance contributions to the residual and Jacobian.
CRSMatrixView< real64, globalIndex const > const m_localMatrix
View on the local CRS matrix.
arrayView1d< integer const > const m_elemGhostRank
View on the ghost ranks.
localIndex const m_iwelemControl
Index of the element where the control is enforced.
GEOS_HOST_DEVICE void setup(localIndex const ei, StackVariables &stack) const
Performs the setup phase for the kernel.
static void createAndLaunch(integer const numComps, real64 const dt, globalIndex const rankOffset, BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags, string const dofKey, WellControls const &wellControls, WellElementSubRegion const &subRegion, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs)
Create a new kernel and launch.
Define the interface for the assembly kernel in charge of flux terms.
arrayView1d< globalIndex const > const m_wellElemDofNumber
Reference to the degree-of-freedom numbers.
arrayView1d< real64 const > const m_injection
Injection stream composition.
GEOS_HOST_DEVICE void setup(localIndex const iconn, StackVariables &stack) const
Performs the setup phase for the kernel.
GEOS_HOST_DEVICE void computeFlux(localIndex const iwelem, StackVariables &stack, FUNC &&compFluxKernelOp=NoOpFunc{}) const
Compute the local flux contributions to the residual and Jacobian.
GEOS_HOST_DEVICE void complete(localIndex const iconn, StackVariables &stack) const
Performs the setup phase for the kernel.
integer const m_rankOffset
Rank offset for calculating row/col Jacobian indices.
arrayView1d< localIndex const > const m_nextWellElemIndex
Next element index, needed since iterating over element nodes, not edges.
arrayView2d< real64 const, compflow::USD_COMP > const m_wellElemCompFrac
Element component fraction.
static void launch(localIndex const numElements, KERNEL_TYPE const &kernelComponent)
Performs the kernel launch.
FaceBasedAssemblyKernel(real64 const dt, globalIndex const rankOffset, string const wellDofKey, WellControls const &wellControls, WellElementSubRegion const &subRegion, CRSMatrixView< real64, globalIndex const > const &localMatrix, arrayView1d< real64 > const &localRhs, BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags)
Constructor for the kernel interface.
arrayView1d< integer const > const m_elemStatus
View on the well status.
CRSMatrixView< real64, globalIndex const > const m_localMatrix
View on the local CRS matrix.
arrayView3d< real64 const, compflow::USD_COMP_DC > const m_dWellElemCompFrac_dCompDens
Element component fraction derivatives.
static void createAndLaunch(integer const numComp, integer const numDof, integer const targetPhaseIndex, globalIndex const rankOffset, string const &dofKey, arrayView1d< real64 const > const &localResidual, WellElementSubRegion const &subRegion, constitutive::MultiFluidBase const &fluid, WellControls const &wellControls, real64 const time, real64 const dt, real64 const minNormalizer, real64(&residualNorm)[1])
Create a new kernel and launch.
virtual GEOS_HOST_DEVICE void computeLinf(localIndex const iwelem, LinfStackVariables &stack) const override
Compute the local values for the Linf norm.
bool const m_isLocallyOwned
Flag indicating whether the well is locally owned or not.
arrayView3d< real64 const, constitutive::multifluid::USD_PHASE > const m_phaseDens_n
View on phase/total density at the previous converged time step.
localIndex const m_iwelemControl
Index of the element where the control is enforced.
virtual GEOS_HOST_DEVICE void computeL2(localIndex const iwelem, L2StackVariables &stack) const override
Compute the local values and normalizer for the L2 norm.
bool const m_isProducer
Flag indicating whether the well is a producer or an injector.
static isothermalCompositionalMultiphaseBaseKernels::SolutionScalingKernel::StackVariables createAndLaunch(real64 const maxRelativePresChange, real64 const maxAbsolutePresChange, real64 const maxCompFracChange, real64 const maxRelativeCompDensChange, globalIndex const rankOffset, integer const numComp, string const dofKey, ElementSubRegionBase &subRegion, arrayView1d< real64 const > const localSolution)
Create a new kernel and launch.
static void createAndLaunch(integer const numComp, integer const numPhase, ObjectManagerBase &subRegion, constitutive::MultiFluidBase const &fluid)
Create a new kernel and launch.
Define the interface for the property kernel in charge of computing the total mass density.
GEOS_HOST_DEVICE void compute(localIndex const ei, FUNC &&totalMassDensityKernelOp=NoOpFunc{}) const
Compute the total mass density in an element.
static constexpr integer numPhase
Compile time value for the number of phases.
arrayView3d< real64 const, constitutive::multifluid::USD_PHASE > m_phaseMassDens
Views on phase mass densities.
arrayView2d< real64 const, compflow::USD_PHASE > m_phaseVolFrac
Views on phase volume fractions.
static constexpr integer numComp
Compile time value for the number of components.
TotalMassDensityKernel(ObjectManagerBase &subRegion, constitutive::MultiFluidBase const &fluid)
Constructor.
GEOS_DECLTYPE_AUTO_RETURN getReference(LOOKUP_TYPE const &lookup) const
Look up a wrapper and get reference to wrapped object.
Definition: Group.hpp:1275
localIndex size() const
Get the "size" of the group, which determines the number of elements in resizable wrappers.
Definition: Group.hpp:1317
Define the base interface for the property update kernels.
static constexpr integer numComp
Compile time value for the number of components.
Define the base interface for the residual calculations.
real64 const m_minNormalizer
Value used to make sure that normalizers are never zero.
arrayView1d< globalIndex const > const m_dofNumber
View on the dof numbers.
GEOS_HOST_DEVICE integer ghostRank(localIndex const i) const
Getter for the ghost rank.
arrayView1d< real64 const > const m_localResidual
View on the local residual.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
StackArray< T, 2, MAXSIZE > stackArray2d
Alias for 2D stack array.
Definition: DataTypes.hpp:203
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:87
StackArray< T, 1, MAXSIZE > stackArray1d
Alias for 1D stack array.
Definition: DataTypes.hpp:187
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
Definition: DataTypes.hpp:199
ArraySlice< T, 3, USD > arraySlice3d
Alias for 3D array slice.
Definition: DataTypes.hpp:215
ArrayView< T, 5, USD > arrayView5d
Alias for 5D array view.
Definition: DataTypes.hpp:243
double real64
64-bit floating point type.
Definition: DataTypes.hpp:98
GEOS_LOCALINDEX_TYPE localIndex
Local index type (for indexing objects within an MPI partition).
Definition: DataTypes.hpp:84
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
Definition: DataTypes.hpp:183
LvArray::CRSMatrixView< T, COL_INDEX, INDEX_TYPE const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:309
ArrayView< T, 4, USD > arrayView4d
Alias for 4D array view.
Definition: DataTypes.hpp:227
ArrayView< T, 2, USD > arrayView2d
Alias for 2D array view.
Definition: DataTypes.hpp:195
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Array< T, 1 > array1d
Alias for 1D array.
Definition: DataTypes.hpp:175
ArrayView< T, 3, USD > arrayView3d
Alias for 3D array view.
Definition: DataTypes.hpp:211
Kernel variables (dof numbers, jacobian and residual) located on the stack.
Kernel variables (dof numbers, jacobian and residual) located on the stack.
real64 localJacobian[numEqn][numDof]
C-array storage for the element local Jacobian matrix (all equations except constraint and momentum)
globalIndex dofIndices[numDof]
Indices of the matrix rows/columns corresponding to the dofs in this element.
real64 localResidual[numEqn]
C-array storage for the element local residual vector (all equations except constraint and momentum)
stackArray2d< real64, maxNumElems *numEqn *maxStencilSize *CP_Deriv::nDer > localFluxJacobian
Storage for the face local Jacobian matrix dC dP dT.
GEOS_HOST_DEVICE StackVariables(localIndex const size)
Constructor for the stack variables.
stackArray2d< real64, maxNumElems *numEqn *maxStencilSize > localFluxJacobian_dQ
Storage for the face local Jacobian matrix dQ only.
stackArray1d< globalIndex, maxNumElems *numDof > dofColIndices
Indices of the matrix rows/columns corresponding to the dofs in this face.
stackArray1d< real64, maxNumElems *numEqn > localFlux
Storage for the face local residual vector (all mass bal equations)
ElementRegionManager::ElementViewConst< VIEWTYPE > ElementViewConst
The type for element-based non-constitutive data parameters. Consists entirely of ArrayView's.