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"
40 #include "physicsSolvers/fluidFlow/wells/WellControls.hpp"
42 
43 namespace geos
44 {
45 
46 
47 namespace compositionalMultiphaseWellKernels
48 {
49 
50 static constexpr real64 minDensForDivision = 1e-10;
51 
52 // tag to access well and reservoir elements in perforation rates computation
54 {
55  static constexpr integer RES = 0;
56  static constexpr integer WELL = 1;
57 };
58 
59 // tag to access the next and current well elements of a connection
60 struct ElemTag
61 {
62  static constexpr integer CURRENT = 0;
63  static constexpr integer NEXT = 1;
64 };
65 
66 // define the column offset of the derivatives
67 struct ColOffset
68 {
69  static constexpr integer DPRES = 0;
70  static constexpr integer DCOMP = 1;
71 };
72 
73 template< integer NC, integer IS_THERMAL >
75 
76 template< integer NC >
77 struct ColOffset_WellJac< NC, 0 >
78 {
79  static constexpr integer dP = 0;
80  static constexpr integer dC = 1;
81  static constexpr integer dQ = dC + NC;
82  static integer constexpr nDer = dQ + 1;
83 
84 };
85 
86 template< integer NC >
87 struct ColOffset_WellJac< NC, 1 >
88 {
89  static constexpr integer dP = 0;
90  static constexpr integer dC = 1;
91  static constexpr integer dQ = dC + NC;
92  static constexpr integer dT = dQ+1;
94  static integer constexpr nDer = dT + 1;
95 };
96 
97 // define the row offset of the residual equations
98 struct RowOffset
99 {
100  static constexpr integer CONTROL = 0;
101  static constexpr integer MASSBAL = 1;
102 };
103 
104 template< integer NC, integer IS_THERMAL >
106 
107 template< integer NC >
108 struct RowOffset_WellJac< NC, 0 >
109 {
110  static constexpr integer CONTROL = 0;
111  static constexpr integer MASSBAL = 1;
112  static constexpr integer VOLBAL = MASSBAL + NC;
113  static constexpr integer nEqn = VOLBAL+1;
114 };
115 
116 template< integer NC >
117 struct RowOffset_WellJac< NC, 1 >
118 {
119  static constexpr integer CONTROL = 0;
120  static constexpr integer MASSBAL = 1;
121  static constexpr integer VOLBAL = MASSBAL + NC;
122  static constexpr integer ENERGYBAL = VOLBAL+1;
123  static constexpr integer nEqn = ENERGYBAL+1;
124 
125 };
126 /******************************** ControlEquationHelper ********************************/
128 {
131 
133  inline
134  static
135  void
136  switchControl( bool const isProducer,
137  WellControls::Control const & inputControl,
138  WellControls::Control const & currentControl,
139  integer const phasePhaseIndex,
140  real64 const & targetBHP,
141  real64 const & targetPhaseRate,
142  real64 const & targetTotalRate,
143  real64 const & targetMassRate,
144  real64 const & currentBHP,
145  arrayView1d< real64 const > const & currentPhaseVolRate,
146  real64 const & currentTotalVolRate,
147  real64 const & currentMassRate,
148  WellControls::Control & newControl );
149 
150  template< integer NC, integer IS_THERMAL >
152  inline
153  static void
154  compute( globalIndex const rankOffset,
155  WellControls::Control const currentControl,
156  integer const targetPhaseIndex,
157  real64 const & targetBHP,
158  real64 const & targetPhaseRate,
159  real64 const & targetTotalRate,
160  real64 const & targetMassRate,
161  real64 const & currentBHP,
162  arrayView1d< real64 const > const & dCurrentBHP,
163  arrayView1d< real64 const > const & currentPhaseVolRate,
164  arrayView2d< real64 const > const & dCurrentPhaseVolRate,
165  real64 const & currentTotalVolRate,
166  arrayView1d< real64 const > const & dCurrentTotalVolRate,
167  real64 const & massDensity,
168  globalIndex const dofNumber,
169  CRSMatrixView< real64, globalIndex const > const & localMatrix,
170  arrayView1d< real64 > const & localRhs );
171 
172 };
173 
174 /******************************** PressureRelationKernel ********************************/
175 
177 {
178  using Deriv = constitutive::multifluid::DerivativeOffset;
182 
183  template< integer NC, integer IS_THERMAL >
185  inline
186  static void
187  compute( real64 const & gravCoef,
188  real64 const & gravCoefNext,
189  real64 const & pres,
190  real64 const & presNext,
191  real64 const & totalMassDens,
192  real64 const & totalMassDensNext,
195  real64 & localPresRel,
196  real64 ( &localPresRelJacobian )[2*(NC+1+IS_THERMAL)] );
197 
198  template< integer NC, integer IS_THERMAL >
199  static void
200  launch( localIndex const size,
201  globalIndex const rankOffset,
202  bool const isLocallyOwned,
203  localIndex const iwelemControl,
204  integer const targetPhaseIndex,
205  WellControls const & wellControls,
206  real64 const & timeAtEndOfStep,
207  arrayView1d< globalIndex const > const & wellElemDofNumber,
208  arrayView1d< real64 const > const & wellElemGravCoef,
209  arrayView1d< localIndex const > const & nextWellElemIndex,
210  arrayView1d< real64 const > const & wellElemPressure,
211  arrayView1d< real64 const > const & wellElemTotalMassDens,
212  arrayView2d< real64 const, compflow::USD_FLUID_DC > const & dWellElemTotalMassDens,
213  bool & controlHasSwitched,
214  CRSMatrixView< real64, globalIndex const > const & localMatrix,
215  arrayView1d< real64 > const & localRhs );
216 
217 };
218 
219 /******************************** VolumeBalanceKernel ********************************/
220 
222 {
223 
226 
227  template< integer NC >
229  inline
230  static void
231  compute( integer const numPhases,
232  real64 const & volume,
235  real64 & localVolBalance,
236  real64 ( &localVolBalanceJacobian )[NC+1] );
237 
238  template< integer NC >
239  static void
240  launch( localIndex const size,
241  integer const numPhases,
242  globalIndex const rankOffset,
243  arrayView1d< globalIndex const > const & wellElemDofNumber,
244  arrayView1d< integer const > const & wellElemGhostRank,
245  arrayView2d< real64 const, compflow::USD_PHASE > const & wellElemPhaseVolFrac,
246  arrayView3d< real64 const, compflow::USD_PHASE_DC > const & dWellElemPhaseVolFrac,
247  arrayView1d< real64 const > const & wellElemVolume,
248  CRSMatrixView< real64, globalIndex const > const & localMatrix,
249  arrayView1d< real64 > const & localRhs );
250 
251 };
252 
253 /******************************** PresTempCompFracInitializationKernel ********************************/
254 
256 {
257 
258  using CompFlowAccessors =
259  StencilAccessors< fields::flow::pressure,
260  fields::flow::temperature,
261  fields::flow::globalCompDensity,
262  fields::flow::phaseVolumeFraction >;
263 
264  using MultiFluidAccessors =
265  StencilMaterialAccessors< constitutive::MultiFluidBase,
266  fields::multifluid::phaseMassDensity >;
267 
268 
276  template< typename VIEWTYPE >
278 
279  static void
280  launch( localIndex const perforationSize,
281  localIndex const subRegionSize,
282  integer const numComponents,
283  integer const numPhases,
284  localIndex const numPerforations,
285  WellControls const & wellControls,
286  real64 const & currentTime,
292  arrayView1d< localIndex const > const & resElementRegion,
293  arrayView1d< localIndex const > const & resElementSubRegion,
294  arrayView1d< localIndex const > const & resElementIndex,
295  arrayView1d< real64 const > const & perfGravCoef,
296  arrayView1d< real64 const > const & wellElemGravCoef,
297  arrayView1d< real64 > const & wellElemPres,
298  arrayView1d< real64 > const & wellElemTemp,
299  arrayView2d< real64, compflow::USD_COMP > const & wellElemCompFrac );
300 
301 };
302 
303 /******************************** CompDensInitializationKernel ********************************/
304 
306 {
307 
308  static void
309  launch( localIndex const subRegionSize,
310  integer const numComponents,
311  arrayView2d< real64 const, compflow::USD_COMP > const & wellElemCompFrac,
313  arrayView2d< real64, compflow::USD_COMP > const & wellElemCompDens );
314 
315 };
316 
317 /******************************** RateInitializationKernel ********************************/
318 
320 {
321 
322  static void
323  launch( localIndex const subRegionSize,
324  integer const targetPhaseIndex,
325  WellControls const & wellControls,
326  real64 const & currentTime,
329  arrayView1d< real64 > const & connRate );
330 
331 };
332 
333 
334 /******************************** TotalMassDensityKernel ****************************/
335 
342 template< integer NUM_COMP, integer NUM_PHASE >
344 {
345 public:
346 
348  using Base::numComp;
349 
351  static constexpr integer numPhase = NUM_PHASE;
352 
359  constitutive::MultiFluidBase const & fluid )
360  : Base(),
361  m_phaseVolFrac( subRegion.getField< fields::well::phaseVolumeFraction >() ),
362  m_dPhaseVolFrac( subRegion.getField< fields::well::dPhaseVolumeFraction >() ),
363  m_dCompFrac_dCompDens( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
364  m_phaseMassDens( fluid.phaseMassDensity() ),
365  m_dPhaseMassDens( fluid.dPhaseMassDensity() ),
366  m_totalMassDens( subRegion.getField< fields::well::totalMassDensity >() ),
367  m_dTotalMassDens( subRegion.getField< fields::well::dTotalMassDensity >() )
368  {}
369 
376  template< typename FUNC = NoOpFunc >
378  inline
379  void compute( localIndex const ei,
380  FUNC && totalMassDensityKernelOp = NoOpFunc{} ) const
381  {
382  using Deriv = constitutive::multifluid::DerivativeOffset;
383 
384  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
385  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
386  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
387  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseMassDens = m_phaseMassDens[ei][0];
388  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseMassDens = m_dPhaseMassDens[ei][0];
389  real64 & totalMassDens = m_totalMassDens[ei];
390  arraySlice1d< real64, compflow::USD_FLUID_DC - 1 > dTotalMassDens = m_dTotalMassDens[ei];
391 
392  real64 dMassDens_dC[numComp]{};
393 
394  totalMassDens = 0.0;
395 
396  dTotalMassDens[Deriv::dP]=0.0;
397  for( integer ic = 0; ic < numComp; ++ic )
398  {
399  dTotalMassDens[Deriv::dC+ic]=0.0;
400  }
401 
402  for( integer ip = 0; ip < numPhase; ++ip )
403  {
404  totalMassDens += phaseVolFrac[ip] * phaseMassDens[ip];
405  dTotalMassDens[Deriv::dP] += dPhaseVolFrac[ip][Deriv::dP] * phaseMassDens[ip] + phaseVolFrac[ip] * dPhaseMassDens[ip][Deriv::dP];
406 
407  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseMassDens[ip], dMassDens_dC, Deriv::dC );
408  for( integer ic = 0; ic < numComp; ++ic )
409  {
410  dTotalMassDens[Deriv::dC+ic] += dPhaseVolFrac[ip][Deriv::dC+ic] * phaseMassDens[ip]
411  + phaseVolFrac[ip] * dMassDens_dC[ic];
412  }
413 
414  totalMassDensityKernelOp( ip ); //, phaseVolFrac, dTotalMassDens_dPres, dTotalMassDens_dCompDens );
415  }
416 
417  }
418 
419 protected:
420 
421  // inputs
422 
427 
431 
432  // outputs
433 
437 
438 
439 };
440 
445 {
446 public:
447 
456  template< typename POLICY >
457  static void
458  createAndLaunch( integer const numComp,
459  integer const numPhase,
460  ObjectManagerBase & subRegion,
461  constitutive::MultiFluidBase const & fluid )
462  {
463  if( numPhase == 2 )
464  {
465  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
466  {
467  integer constexpr NUM_COMP = NC();
468  TotalMassDensityKernel< NUM_COMP, 2 > kernel( subRegion, fluid );
469  TotalMassDensityKernel< NUM_COMP, 2 >::template launch< POLICY >( subRegion.size(), kernel );
470  } );
471  }
472  else if( numPhase == 3 )
473  {
474  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
475  {
476  integer constexpr NUM_COMP = NC();
477  TotalMassDensityKernel< NUM_COMP, 3 > kernel( subRegion, fluid );
478  TotalMassDensityKernel< NUM_COMP, 3 >::template launch< POLICY >( subRegion.size(), kernel );
479  } );
480  }
481  }
482 };
483 
484 
485 /******************************** ResidualNormKernel ********************************/
486 
491 {
492 public:
493 
495  using Base::m_minNormalizer;
496  using Base::m_rankOffset;
497  using Base::m_localResidual;
498  using Base::m_dofNumber;
499 
500  ResidualNormKernel( globalIndex const rankOffset,
501  arrayView1d< real64 const > const & localResidual,
502  arrayView1d< globalIndex const > const & dofNumber,
504  integer const numComp,
505  integer const numDof,
506  integer const targetPhaseIndex,
507  WellElementSubRegion const & subRegion,
508  constitutive::MultiFluidBase const & fluid,
509  WellControls const & wellControls,
510  real64 const timeAtEndOfStep,
511  real64 const dt,
512  real64 const minNormalizer )
513  : Base( rankOffset,
514  localResidual,
515  dofNumber,
516  ghostRank,
517  minNormalizer ),
518  m_numComp( numComp ),
519  m_numDof( numDof ),
520  m_targetPhaseIndex( targetPhaseIndex ),
521  m_dt( dt ),
522  m_isLocallyOwned( subRegion.isLocallyOwned() ),
524  m_isProducer( wellControls.isProducer() ),
525  m_currentControl( wellControls.getControl() ),
526  m_targetBHP( wellControls.getTargetBHP( timeAtEndOfStep ) ),
527  m_targetTotalRate( wellControls.getTargetTotalRate( timeAtEndOfStep ) ),
528  m_targetPhaseRate( wellControls.getTargetPhaseRate( timeAtEndOfStep ) ),
529  m_targetMassRate( wellControls.getTargetMassRate( timeAtEndOfStep ) ),
530  m_volume( subRegion.getElementVolume() ),
531  m_phaseDens_n( fluid.phaseDensity_n() ),
532  m_totalDens_n( fluid.totalDensity_n() )
533  {}
534 
536  virtual void computeLinf( localIndex const iwelem,
537  LinfStackVariables & stack ) const override
538  {
540 
541  real64 normalizer = 0.0;
542  for( integer idof = 0; idof < m_numDof; ++idof )
543  {
544 
545  // Step 1: compute a normalizer for the control or pressure equation
546 
547  // for the control equation, we distinguish two cases
548  if( idof == ROFFSET::CONTROL )
549  {
550 
551  // for the top well element, normalize using the current control
552  if( m_isLocallyOwned && iwelem == m_iwelemControl )
553  {
555  {
556  // the residual entry is in pressure units
557  normalizer = m_targetBHP;
558  }
560  {
561  // the residual entry is in volume / time units
562  normalizer = LvArray::math::max( LvArray::math::abs( m_targetTotalRate ), m_minNormalizer );
563  }
565  {
566  // the residual entry is in volume / time units
567  normalizer = LvArray::math::max( LvArray::math::abs( m_targetPhaseRate ), m_minNormalizer );
568  }
570  {
571  // the residual entry is in volume / time units
572  normalizer = LvArray::math::max( LvArray::math::abs( m_targetMassRate ), m_minNormalizer );
573  }
574  }
575  // for the pressure difference equation, always normalize by the BHP
576  else
577  {
578  normalizer = m_targetBHP;
579  }
580  }
581  // Step 2: compute a normalizer for the mass balance equations
582  else if( idof >= ROFFSET::MASSBAL && idof < ROFFSET::MASSBAL + m_numComp )
583  {
584  if( m_isProducer ) // only PHASEVOLRATE is supported for now
585  {
586  // the residual is in mass units
587  normalizer = m_dt * LvArray::math::abs( m_targetPhaseRate ) * m_phaseDens_n[iwelem][0][m_targetPhaseIndex];
588  }
589  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
590  {
592  {
593  normalizer = m_dt * LvArray::math::abs( m_targetMassRate );
594  }
595  else
596  {
597  // the residual is in mass units
598  normalizer = m_dt * LvArray::math::abs( m_targetTotalRate ) * m_totalDens_n[iwelem][0];
599  }
600 
601  }
602 
603  // to make sure that everything still works well if the rate is zero, we add this check
604  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] * m_totalDens_n[iwelem][0] );
605  }
606  // Step 3: compute a normalizer for the volume balance equations
607  else if( idof == ROFFSET::MASSBAL + m_numComp )
608  {
609  if( m_isProducer ) // only PHASEVOLRATE is supported for now
610  {
611  // the residual is in volume units
612  normalizer = m_dt * LvArray::math::abs( m_targetPhaseRate );
613  }
614  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
615  {
617  {
618  normalizer = m_dt * LvArray::math::abs( m_targetMassRate/ m_totalDens_n[iwelem][0] );
619  }
620  else
621  {
622  normalizer = m_dt * LvArray::math::abs( m_targetTotalRate );
623  }
624 
625  }
626 
627  }
628 
629  // to make sure that everything still works well if the rate is zero, we add this check
630  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] );
631 
632  // Step 4: compute the contribution to the residual
633  real64 const val = LvArray::math::abs( m_localResidual[stack.localRow + idof] ) / normalizer;
634  if( val > stack.localValue[0] )
635  {
636  stack.localValue[0] = val;
637  }
638  }
639  }
640 
642  virtual void computeL2( localIndex const iwelem,
643  L2StackVariables & stack ) const override
644  {
645  GEOS_UNUSED_VAR( iwelem, stack );
646  GEOS_ERROR( "The L2 norm is not implemented for CompositionalMultiphaseWell" );
647  }
648 
649 
650 protected:
651 
654 
657 
660 
662  real64 const m_dt;
663 
665  bool const m_isLocallyOwned;
666 
669 
671  bool const m_isProducer;
672 
675  real64 const m_targetBHP;
676  real64 const m_targetTotalRate;
677  real64 const m_targetPhaseRate;
678  real64 const m_targetMassRate;
679 
682 
686 
687 };
688 
693 {
694 public:
695 
712  template< typename POLICY >
713  static void
714  createAndLaunch( integer const numComp,
715  integer const numDof,
716  integer const targetPhaseIndex,
717  globalIndex const rankOffset,
718  string const & dofKey,
719  arrayView1d< real64 const > const & localResidual,
720  WellElementSubRegion const & subRegion,
721  constitutive::MultiFluidBase const & fluid,
722  WellControls const & wellControls,
723  real64 const timeAtEndOfStep,
724  real64 const dt,
725  real64 const minNormalizer,
726  real64 (& residualNorm)[1] )
727  {
728  arrayView1d< globalIndex const > const dofNumber = subRegion.getReference< array1d< globalIndex > >( dofKey );
729  arrayView1d< integer const > const ghostRank = subRegion.ghostRank();
730 
731  ResidualNormKernel kernel( rankOffset, localResidual, dofNumber, ghostRank,
732  numComp, numDof, targetPhaseIndex, subRegion, fluid, wellControls, timeAtEndOfStep, dt, minNormalizer );
733  ResidualNormKernel::launchLinf< POLICY >( subRegion.size(), kernel, residualNorm );
734  }
735 
736 };
737 
738 /******************************** SolutionScalingKernel ********************************/
739 
744 {
745 public:
746 
760  template< typename POLICY >
762  createAndLaunch( real64 const maxRelativePresChange,
763  real64 const maxAbsolutePresChange,
764  real64 const maxCompFracChange,
765  real64 const maxRelativeCompDensChange,
766  globalIndex const rankOffset,
767  integer const numComp,
768  string const dofKey,
769  ElementSubRegionBase & subRegion,
770  arrayView1d< real64 const > const localSolution )
771  {
772  arrayView1d< real64 const > const pressure =
773  subRegion.getField< fields::well::pressure >();
775  subRegion.getField< fields::well::globalCompDensity >();
776  arrayView1d< real64 > pressureScalingFactor =
777  subRegion.getField< fields::well::pressureScalingFactor >();
778  arrayView1d< real64 > compDensScalingFactor =
779  subRegion.getField< fields::well::globalCompDensityScalingFactor >();
780  isothermalCompositionalMultiphaseBaseKernels::
781  SolutionScalingKernel kernel( maxRelativePresChange, maxAbsolutePresChange, maxCompFracChange, maxRelativeCompDensChange, rankOffset,
782  numComp, dofKey, subRegion, localSolution, pressure, compDens, pressureScalingFactor, compDensScalingFactor );
783  return isothermalCompositionalMultiphaseBaseKernels::
784  SolutionScalingKernel::
785  launch< POLICY >( subRegion.size(), kernel );
786  }
787 
788 };
789 
790 /******************************** SolutionCheckKernel ********************************/
791 
796 {
797 public:
798 
810  template< typename POLICY >
812  createAndLaunch( integer const allowCompDensChopping,
813  CompositionalMultiphaseFVM::ScalingType const scalingType,
814  real64 const scalingFactor,
815  arrayView1d< real64 const > const pressure,
817  arrayView1d< real64 > pressureScalingFactor,
818  arrayView1d< real64 > compDensScalingFactor,
819  globalIndex const rankOffset,
820  integer const numComp,
821  string const dofKey,
822  ElementSubRegionBase & subRegion,
823  arrayView1d< real64 const > const localSolution )
824  {
825 
826  isothermalCompositionalMultiphaseBaseKernels::
827  SolutionCheckKernel kernel( allowCompDensChopping, 0, scalingType, scalingFactor,
828  pressure, compDens, pressureScalingFactor, compDensScalingFactor, rankOffset,
829  numComp, dofKey, subRegion, localSolution );
830  return isothermalCompositionalMultiphaseBaseKernels::
831  SolutionCheckKernel::
832  launch< POLICY >( subRegion.size(), kernel );
833  }
834 
835 };
836 
837 /******************************** ElementBasedAssemblyKernel ********************************/
838 
845 template< integer NUM_COMP, integer IS_THERMAL >
847 {
848 public:
851 
852  // Well jacobian column and row indicies
853  // tjb - change NUM_DOF to IS_THERMAL
854  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NUM_COMP, IS_THERMAL >;
858  static constexpr integer numComp = NUM_COMP;
859 
861  static constexpr integer numDof = NUM_COMP + 1 + IS_THERMAL;
862 
864  static constexpr integer numEqn = NUM_COMP + 1 + IS_THERMAL;
865 
866 
879  integer const isProducer,
880  globalIndex const rankOffset,
881  string const dofKey,
882  WellElementSubRegion const & subRegion,
883  constitutive::MultiFluidBase const & fluid,
884  CRSMatrixView< real64, globalIndex const > const & localMatrix,
885  arrayView1d< real64 > const & localRhs,
886  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const kernelFlags )
887  : m_numPhases( numPhases ),
888  m_isProducer( isProducer ),
889  m_rankOffset( rankOffset ),
890  m_iwelemControl( subRegion.getTopWellElementIndex() ),
891  m_dofNumber( subRegion.getReference< array1d< globalIndex > >( dofKey ) ),
892  m_elemGhostRank( subRegion.ghostRank() ),
893  m_volume( subRegion.getElementVolume() ),
894  m_dCompFrac_dCompDens( subRegion.getField< fields::flow::dGlobalCompFraction_dGlobalCompDensity >() ),
895  m_phaseVolFrac_n( subRegion.getField< fields::flow::phaseVolumeFraction_n >() ),
896  m_phaseVolFrac( subRegion.getField< fields::flow::phaseVolumeFraction >() ),
897  m_dPhaseVolFrac( subRegion.getField< fields::flow::dPhaseVolumeFraction >() ),
898  m_phaseDens_n( fluid.phaseDensity_n() ),
899  m_phaseDens( fluid.phaseDensity() ),
900  m_dPhaseDens( fluid.dPhaseDensity() ),
901  m_phaseCompFrac_n( fluid.phaseCompFraction_n() ),
902  m_phaseCompFrac( fluid.phaseCompFraction() ),
903  m_dPhaseCompFrac( fluid.dPhaseCompFraction() ),
904  m_compDens( subRegion.getField< fields::flow::globalCompDensity >() ),
905  m_compDens_n( subRegion.getField< fields::flow::globalCompDensity_n >() ),
906  m_localMatrix( localMatrix ),
907  m_localRhs( localRhs ),
908  m_kernelFlags( kernelFlags )
909  {}
910 
916  {
917 public:
918 
919  // volume information (used by both accumulation and volume balance)
920  real64 volume = 0.0;
921 
922  // Residual information
923 
926 
928  globalIndex dofIndices[numDof]{}; // NC compdens + P + thermal
929  globalIndex eqnRowIndices[numDof]{};
930  globalIndex dofColIndices[numDof]{};
931 
934 
937 
938  };
939 
946  integer elemGhostRank( localIndex const ei ) const
947  { return m_elemGhostRank( ei ); }
948 
949 
956  void setup( localIndex const ei,
957  StackVariables & stack ) const
958  {
959  // initialize the volume
960  stack.volume = m_volume[ei];
961 
962  // Note row/col indices needed to be consistent with layout of stack.localJacobian
963  // Setup row equation indices for this element ( mass + vol + thermal if valid)
964 
965  // 1) Mass Balance
966  for( integer ic = 0; ic < numComp; ++ic )
967  {
968  stack.eqnRowIndices[ic] = m_dofNumber[ei] + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
969  }
970 // 2) Volume Balance
971  stack.eqnRowIndices[numComp] = m_dofNumber[ei] + WJ_ROFFSET::VOLBAL - m_rankOffset;
972  // 3) Energy Balance
973  if constexpr ( IS_THERMAL )
974  {
975  stack.eqnRowIndices[numComp+1] = m_dofNumber[ei] + WJ_ROFFSET::ENERGYBAL - m_rankOffset;
976  }
977  // Setup equation column indices for this element ( P + COMPDENS + THERMAL if valid)
978  stack.dofColIndices[0] = m_dofNumber[ei] + WJ_COFFSET::dP;
979  for( integer ic = 0; ic < numComp; ++ic )
980  {
981  stack.dofColIndices[ic+1] = m_dofNumber[ei] + WJ_COFFSET::dC+ic;
982  }
983  if constexpr ( IS_THERMAL )
984  {
985  stack.dofColIndices[numComp+1] = m_dofNumber[ei] + WJ_COFFSET::dT;
986  }
987  if( 1 )
988  for( integer jc = 0; jc < numEqn; ++jc )
989  {
990  stack.localResidual[jc] = 0.0;
991  for( integer ic = 0; ic < numDof; ++ic )
992  {
993  stack.localJacobian[jc][ic] = 0.0;
994  }
995 
996  }
997 
998  }
999 
1007  template< typename FUNC = NoOpFunc >
1010  StackVariables & stack,
1011  FUNC && phaseAmountKernelOp = NoOpFunc{} ) const
1012  {
1013 
1014  using Deriv = constitutive::multifluid::DerivativeOffset;
1015 
1016  // construct the slices for variables accessed multiple times
1017  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
1018 
1019  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac_n = m_phaseVolFrac_n[ei];
1020  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
1021  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
1022 
1023  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens_n = m_phaseDens_n[ei][0];
1024  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens = m_phaseDens[ei][0];
1025  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseDens = m_dPhaseDens[ei][0];
1026 
1027  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac_n = m_phaseCompFrac_n[ei][0];
1028  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac = m_phaseCompFrac[ei][0];
1029  arraySlice3d< real64 const, constitutive::multifluid::USD_PHASE_COMP_DC - 2 > dPhaseCompFrac = m_dPhaseCompFrac[ei][0];
1030 
1031  // temporary work arrays
1032  real64 dPhaseAmount[FLUID_PROP_COFFSET::nDer]{};
1033  real64 dPhaseAmount_dC[numComp]{};
1034  real64 dPhaseCompFrac_dC[numComp]{};
1035 
1036  // sum contributions to component accumulation from each phase
1037  for( integer ip = 0; ip < m_numPhases; ++ip )
1038  {
1039  real64 const phaseAmount = stack.volume * phaseVolFrac[ip] * phaseDens[ip];
1040  real64 const phaseAmount_n = stack.volume * phaseVolFrac_n[ip] * phaseDens_n[ip];
1041  //remove tjb
1042  real64 const dPhaseAmount_dP = stack.volume * ( dPhaseVolFrac[ip][Deriv::dP] * phaseDens[ip]
1043  + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dP] );
1044  dPhaseAmount[FLUID_PROP_COFFSET::dP]=stack.volume * ( dPhaseVolFrac[ip][Deriv::dP] * phaseDens[ip]
1045  + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dP] );
1046 
1047  // assemble density dependence
1048  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], dPhaseAmount_dC, Deriv::dC );
1049  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], &dPhaseAmount[FLUID_PROP_COFFSET::dC], Deriv::dC );
1050  for( integer jc = 0; jc < numComp; ++jc )
1051  {
1052  dPhaseAmount_dC[jc] = dPhaseAmount_dC[jc] * phaseVolFrac[ip]
1053  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1054  dPhaseAmount_dC[jc] *= stack.volume;
1055  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] = dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] * phaseVolFrac[ip]
1056  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1057  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] *= stack.volume;
1058  }
1059 // tjb- remove when safe
1060  for( integer ic = 0; ic < numComp; ic++ )
1061  {
1062  assert( fabs( dPhaseAmount[FLUID_PROP_COFFSET::dC+ic] -dPhaseAmount_dC[ic] ) < FLT_EPSILON );
1063 
1064  }
1065  // ic - index of component whose conservation equation is assembled
1066  // (i.e. row number in local matrix)
1067  for( integer ic = 0; ic < numComp; ++ic )
1068  {
1069  real64 const phaseCompAmount = phaseAmount * phaseCompFrac[ip][ic];
1070  real64 const phaseCompAmount_n = phaseAmount_n * phaseCompFrac_n[ip][ic];
1071 
1072  real64 const dPhaseCompAmount_dP = dPhaseAmount_dP * phaseCompFrac[ip][ic]
1073  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dP];
1074 
1075  stack.localResidual[ic] += phaseCompAmount - phaseCompAmount_n;
1076  stack.localJacobian[ic][0] += dPhaseCompAmount_dP;
1077 
1078  // jc - index of component w.r.t. whose compositional var the derivative is being taken
1079  // (i.e. col number in local matrix)
1080 
1081  // assemble phase composition dependence
1082  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseCompFrac[ip][ic], dPhaseCompFrac_dC, Deriv::dC );
1083  for( integer jc = 0; jc < numComp; ++jc )
1084  {
1085  real64 const dPhaseCompAmount_dC = dPhaseCompFrac_dC[jc] * phaseAmount
1086  + phaseCompFrac[ip][ic] * dPhaseAmount[FLUID_PROP_COFFSET::dC+jc];
1087 
1088  stack.localJacobian[ic][jc + 1] += dPhaseCompAmount_dC;
1089  }
1090  }
1091  if constexpr ( IS_THERMAL )
1092  {
1093  dPhaseAmount[FLUID_PROP_COFFSET::dT] = stack.volume * (dPhaseVolFrac[ip][Deriv::dT] * phaseDens[ip] + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dT] );
1094  for( integer ic = 0; ic < numComp; ++ic )
1095  {
1096  // assemble the derivatives of the component mass balance equations with respect to temperature
1097  stack.localJacobian[ic][numComp+1] += dPhaseAmount[FLUID_PROP_COFFSET::dT] * phaseCompFrac[ip][ic]
1098  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dT];
1099  }
1100  }
1101  // call the lambda in the phase loop to allow the reuse of the phase amounts and their derivatives
1102  // possible use: assemble accumulation term of the energy equation for this phase
1103  phaseAmountKernelOp( ip, phaseAmount, phaseAmount_n, dPhaseAmount );
1104 
1105  }
1106 
1107  // check zero diagonal (works only in debug)
1108  /*
1109  for( integer ic = 0; ic < numComp; ++ic )
1110  {
1111  GEOS_ASSERT_MSG ( LvArray::math::abs( stack.localJacobian[ic][ic] ) > minDensForDivision,
1112  GEOS_FMT( "Zero diagonal in Jacobian: equation {}, value = {}", ic, stack.localJacobian[ic][ic] ) );
1113  }
1114  */
1115  }
1116 
1117 
1128  StackVariables & stack ) const
1129  {
1130  using Deriv = constitutive::multifluid::DerivativeOffset;
1131 
1132  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
1133  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
1134 
1135  real64 oneMinusPhaseVolFracSum = 1.0;
1136 
1137  // sum contributions to component accumulation from each phase
1138 // Note localJacobian stores equation balances in order of component/vol/enerqy
1139  // These are mapped to solver orderings with indicies setup in stack variables
1140  for( integer ip = 0; ip < m_numPhases; ++ip )
1141  {
1142  oneMinusPhaseVolFracSum -= phaseVolFrac[ip];
1143  stack.localJacobian[numComp][0] -= dPhaseVolFrac[ip][Deriv::dP];
1144 
1145  for( integer jc = 0; jc < numComp; ++jc )
1146  {
1147  stack.localJacobian[numComp][jc+1] -= dPhaseVolFrac[ip][Deriv::dC+jc];
1148  }
1149 
1150  if constexpr ( IS_THERMAL)
1151  {
1152  stack.localJacobian[numComp][numComp+1] -= dPhaseVolFrac[ip][Deriv::dT];
1153  }
1154 
1155  }
1156  // scale saturation-based volume balance by pore volume (for better scaling w.r.t. other equations)
1157  stack.localResidual[numComp] = stack.volume * oneMinusPhaseVolFracSum;
1158  for( integer idof = 0; idof < numComp+1+IS_THERMAL; ++idof )
1159  {
1160  stack.localJacobian[numComp][idof] *= stack.volume;
1161  }
1162 
1163  }
1164 
1171  void complete( localIndex const ei, //GEOS_UNUSED_PARAM( ei ),
1172  StackVariables & stack ) const
1173  {
1174  using namespace compositionalMultiphaseUtilities;
1175 
1176  integer const numRows = numComp+1+ IS_THERMAL;
1177 
1178  if constexpr ( IS_THERMAL)
1179  {
1180  if( ei == m_iwelemControl && !m_isProducer )
1181  {
1182  // For top segment energy balance eqn replaced with T(n+1) - T = 0
1183  // No other energy balance derivatives
1184  // Assumption is global index == 0 is top segment with fixed temp BC
1185 
1186  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1187  {
1188  stack.localJacobian[numRows-1][i] = 0.0;
1189  }
1190  // constant Temperature
1191  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1192  stack.localJacobian[i][numRows-1] = 0.0;
1193  stack.localJacobian[numRows-1][numRows-1] = 1.0;
1194 
1195  stack.localResidual[numRows-1]=0.0;
1196  }
1197  }
1198 
1199  if( m_kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) )
1200  {
1201  // apply equation/variable change transformation to the component mass balance equations
1202  real64 work[numComp + 1 + IS_THERMAL]{};
1203  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numComp, numComp+1+ IS_THERMAL, stack.localJacobian, work );
1204  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numComp, stack.localResidual );
1205  }
1206 
1207  // add contribution to residual and jacobian into:
1208  // - the component mass balance equations (i = 0 to i = numComp-1)
1209  // - the volume balance equations (i = numComp)
1210  // note that numDof includes derivatives wrt temperature if this class is derived in ThermalKernels
1211 
1212  for( integer i = 0; i < numRows; ++i )
1213  {
1214  m_localRhs[stack.eqnRowIndices[i]] += stack.localResidual[i];
1215  m_localMatrix.addToRow< serialAtomic >( stack.eqnRowIndices[i],
1216  stack.dofColIndices,
1217  stack.localJacobian[i],
1218  numComp+1+ IS_THERMAL );
1219  }
1220 
1221  }
1222 
1230  template< typename POLICY, typename KERNEL_TYPE >
1231  static void
1232  launch( localIndex const numElems,
1233  KERNEL_TYPE const & kernelComponent )
1234  {
1236 
1237  forAll< POLICY >( numElems, [=] GEOS_HOST_DEVICE ( localIndex const ei )
1238  {
1239  if( kernelComponent.elemGhostRank( ei ) >= 0 )
1240  {
1241  return;
1242  }
1243 
1244  typename KERNEL_TYPE::StackVariables stack;
1245 
1246  kernelComponent.setup( ei, stack );
1247  kernelComponent.computeAccumulation( ei, stack );
1248  kernelComponent.computeVolumeBalance( ei, stack );
1249  kernelComponent.complete( ei, stack );
1250  } );
1251  }
1252 
1253 protected:
1254 
1257 
1260 
1263 
1266 
1269 
1272 
1275 
1278  arrayView2d< real64 const > const m_porosity;
1279  arrayView2d< real64 const > const m_dPoro_dPres;
1280 
1283 
1288 
1293 
1298 
1299  // Views on component densities
1302 
1307 
1308  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const m_kernelFlags;
1309 };
1310 
1311 
1316 {
1317 public:
1330  template< typename POLICY >
1331  static void
1332  createAndLaunch( localIndex const numComps,
1333  localIndex const numPhases,
1334  integer const isProducer,
1335  globalIndex const rankOffset,
1336  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags,
1337  string const dofKey,
1338  WellElementSubRegion const & subRegion,
1339  constitutive::MultiFluidBase const & fluid,
1340  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1341  arrayView1d< real64 > const & localRhs )
1342  {
1343  geos::internal::kernelLaunchSelectorCompThermSwitch( numComps, 0, [&]( auto NC, auto IS_THERMAL )
1344  {
1345  localIndex constexpr NUM_COMP = NC();
1346 
1347  integer constexpr istherm = IS_THERMAL();
1348 
1350  kernel( numPhases, isProducer, rankOffset, dofKey, subRegion, fluid, localMatrix, localRhs, kernelFlags );
1352  launch< POLICY, ElementBasedAssemblyKernel< NUM_COMP, istherm > >( subRegion.size(), kernel );
1353  } );
1354  }
1355 };
1364 template< integer NC, integer IS_THERMAL >
1366 {
1367 public:
1368 
1372 
1373  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1376 
1377  using CP_Deriv = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1379  static constexpr integer numComp = NC;
1380 
1382  static constexpr integer numDof = WJ_COFFSET::nDer;
1383 
1385  static constexpr integer numEqn = NC;
1386 
1387  static constexpr integer maxNumElems = 2;
1388  static constexpr integer maxStencilSize = 2;
1404  globalIndex const rankOffset,
1405  string const wellDofKey,
1406  WellControls const & wellControls,
1407  WellElementSubRegion const & subRegion,
1408  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1409  arrayView1d< real64 > const & localRhs,
1410  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags )
1411  :
1412  m_dt( dt ),
1413  m_rankOffset( rankOffset ),
1414  m_wellElemDofNumber ( subRegion.getReference< array1d< globalIndex > >( wellDofKey ) ),
1415  m_nextWellElemIndex ( subRegion.getReference< array1d< localIndex > >( WellElementSubRegion::viewKeyStruct::nextWellElementIndexString()) ),
1416  m_connRate ( subRegion.getField< fields::well::mixtureConnectionRate >() ),
1417  m_wellElemCompFrac ( subRegion.getField< fields::well::globalCompFraction >() ),
1418  m_dWellElemCompFrac_dCompDens ( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
1419  m_localMatrix( localMatrix ),
1420  m_localRhs ( localRhs ),
1421  m_useTotalMassEquation ( kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) ),
1422  m_isProducer ( wellControls.isProducer() ),
1423  m_injection ( wellControls.getInjectionStream() )
1424  {}
1425 
1427  {
1428 public:
1429 
1437  : stencilSize( size ),
1438  numConnectedElems( 2 ),
1439  dofColIndices( size * numDof )
1440  {}
1441 
1442  // Stencil information
1443  localIndex const stencilSize;
1446 
1447 
1448  // edge indexes
1449  localIndex iwelemUp;
1450  localIndex iwelemNext;
1451  localIndex iwelemCurrent;
1452  globalIndex offsetUp;
1453  globalIndex offsetCurrent;
1454  globalIndex offsetNext;
1455  // Local degrees of freedom and local residual/jacobian
1456 
1459 
1466  };
1467 
1474  inline
1475  void setup( localIndex const iconn,
1476  StackVariables & stack ) const
1477  {
1478  stack.numConnectedElems=2;
1479  if( m_nextWellElemIndex[iconn] <0 )
1480  {
1481  stack.numConnectedElems = 1;
1482  }
1483  stack.localFlux.resize( stack.numConnectedElems*numEqn );
1484  stack.localFluxJacobian.resize( stack.numConnectedElems * numEqn, stack.stencilSize * numDof );
1485  stack.localFluxJacobian_dQ.resize( stack.numConnectedElems * numEqn, 1 );
1486 
1487  }
1488 
1495  inline
1496  void complete( localIndex const iconn,
1497  StackVariables & stack ) const
1498  {
1499  GEOS_UNUSED_VAR( iconn );
1500  using namespace compositionalMultiphaseUtilities;
1501  if( stack.numConnectedElems ==1 )
1502  {
1503  // Setup Jacobian global row indicies
1504  // equations for COMPONENT + ENERGY balances
1505  globalIndex oneSidedEqnRowIndices[numEqn]{};
1506  for( integer ic = 0; ic < NC; ++ic )
1507  {
1508  oneSidedEqnRowIndices[ic] = stack.offsetUp + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1509  }
1510 
1511  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1512  globalIndex oneSidedDofColIndices_dPresCompTempUp[CP_Deriv::nDer]{};
1513  globalIndex oneSidedDofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1514  // Note localFluxJacobian cols are stored using CP_Deriv order (dP dC or dP dT dC)
1515  int ioff=0;
1516  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1517 
1518  if constexpr ( IS_THERMAL )
1519  {
1520  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1521  }
1522  for( integer jdof = 0; jdof < NC; ++jdof )
1523  {
1524  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1525  }
1526  if( m_useTotalMassEquation > 0 )
1527  {
1528  // Apply equation/variable change transformation(s)
1529  real64 work[CP_Deriv::nDer]{};
1530  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, 1, stack.localFluxJacobian_dQ, work );
1531  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, CP_Deriv::nDer, stack.localFluxJacobian, work );
1532  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, stack.localFlux );
1533  }
1534  for( integer i = 0; i < numEqn; ++i )
1535  {
1536  if( oneSidedEqnRowIndices[i] >= 0 && oneSidedEqnRowIndices[i] < m_localMatrix.numRows() )
1537  {
1538  m_localMatrix.addToRow< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1539  &oneSidedDofColIndices_dRate,
1540  stack.localFluxJacobian_dQ[i],
1541  1 );
1542  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1543  oneSidedDofColIndices_dPresCompTempUp,
1544  stack.localFluxJacobian[i],
1545  CP_Deriv::nDer );
1546  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[oneSidedEqnRowIndices[i]], stack.localFlux[i] );
1547  }
1548  }
1549  }
1550  else
1551  {
1552  // Setup Jacobian global row indicies
1553  // equations for COMPONENT + ENERGY balances
1554  globalIndex eqnRowIndices[2*numEqn]{};
1555 
1556  for( integer ic = 0; ic < NC; ++ic )
1557  {
1558  // mass balance equations for all components
1559  eqnRowIndices[TAG::NEXT *numEqn+ic] = stack.offsetNext + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1560  eqnRowIndices[TAG::CURRENT *numEqn+ic] = stack.offsetCurrent + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1561  }
1562 
1563  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1564  globalIndex dofColIndices_dPresCompUp[CP_Deriv::nDer]{};
1565  globalIndex dofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1566 
1567  int ioff=0;
1568  // Indice storage order reflects local jac col storage order CP::Deriv order P T DENS
1569  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1570 
1571  if constexpr ( IS_THERMAL )
1572  {
1573  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1574  }
1575  for( integer jdof = 0; jdof < NC; ++jdof )
1576  {
1577  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1578  }
1579 
1580 
1581  if( m_useTotalMassEquation > 0 )
1582  {
1583  // Apply equation/variable change transformation(s)
1584  real64 work[CP_Deriv::nDer]{};
1585  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, 1, 2, stack.localFluxJacobian_dQ, work );
1586  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, CP_Deriv::nDer, 2, stack.localFluxJacobian, work );
1587  shiftBlockElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, numEqn, 2, stack.localFlux );
1588  }
1589  // Note this updates diag and offdiag
1590  for( integer i = 0; i < 2*NC; ++i )
1591  {
1592  if( eqnRowIndices[i] >= 0 && eqnRowIndices[i] < m_localMatrix.numRows() )
1593  {
1594  m_localMatrix.addToRow< parallelDeviceAtomic >( eqnRowIndices[i],
1595  &dofColIndices_dRate,
1596  stack.localFluxJacobian_dQ[i],
1597  1 );
1598  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( eqnRowIndices[i],
1599  dofColIndices_dPresCompUp,
1600  stack.localFluxJacobian[i],
1601  CP_Deriv::nDer );
1602  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[eqnRowIndices[i]], stack.localFlux[i] );
1603 
1604  }
1605  }
1606  }
1607  }
1608 
1610  inline
1611  void
1612  computeExit( real64 const & dt,
1613  real64 const ( &compFlux )[NC ],
1614  StackVariables & stack,
1615  real64 ( & dCompFlux)[NC][numDof] ) const
1616  {
1617  for( integer ic = 0; ic < NC; ++ic )
1618  {
1619  stack.localFlux[ic] = -dt * compFlux[ic];
1620  // derivative with respect to rate
1621  stack.localFluxJacobian_dQ[ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1622  // derivative with respect to upstream pressure
1623  stack.localFluxJacobian[ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1624  // derivatives with respect to upstream component densities
1625  for( integer jdof = 0; jdof < NC; ++jdof )
1626  {
1627  stack.localFluxJacobian[ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1628  }
1629  if constexpr ( IS_THERMAL )
1630  {
1631  stack.localFluxJacobian[ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1632  }
1633  }
1634  }
1635 
1637  inline
1638  void
1639  compute( real64 const & dt,
1640  real64 const ( &compFlux )[NC ],
1641  StackVariables & stack,
1642  real64 ( & dCompFlux)[(NC )][numDof] ) const
1643  {
1644  // flux terms
1645  for( integer ic = 0; ic < NC; ++ic )
1646  {
1647  stack.localFlux[TAG::NEXT * NC +ic] = dt * compFlux[ic];
1648  stack.localFlux[TAG::CURRENT * NC +ic] = -dt * compFlux[ic];
1649  // derivative with respect to rate
1650  stack.localFluxJacobian_dQ[TAG::NEXT * NC+ ic][0] = dt * dCompFlux[ic][WJ_COFFSET::dQ];
1651  stack.localFluxJacobian_dQ[TAG::CURRENT * NC + +ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1652 
1653  // derivative with respect to upstream pressure
1654  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dP] = dt * dCompFlux[ic][WJ_COFFSET::dP];
1655  stack.localFluxJacobian[TAG::CURRENT * NC+ ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1656 
1657  if constexpr ( IS_THERMAL )
1658  {
1659  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dT] = dt * dCompFlux[ic][WJ_COFFSET::dT];
1660  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1661  }
1662 
1663  // derivatives with respect to upstream component densities
1664  for( integer jdof = 0; jdof < NC; ++jdof )
1665  {
1666  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dC+jdof] = dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1667  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1668  }
1669  }
1670  }
1671 
1679  template< typename FUNC = NoOpFunc >
1681  inline
1682  void computeFlux( localIndex const iwelem,
1683  StackVariables & stack,
1684  FUNC && compFluxKernelOp = NoOpFunc{} ) const
1685  {
1686 
1687  using namespace compositionalMultiphaseUtilities;
1688 
1689  // create local work arrays
1690  real64 compFracUp[NC]{};
1691  real64 compFlux[NC]{};
1692  real64 dComp[NC][NC];
1693  real64 dCompFlux[NC][numDof]{};
1694  for( integer ic = 0; ic < NC; ++ic )
1695  {
1696  for( integer jc = 0; jc < NC; ++jc )
1697  {
1698  dComp[ic][jc]=0.0;
1699  }
1700  }
1701  // Step 1) decide the upwind well element
1702 
1703  /* currentConnRate < 0 flow from iwelem to iwelemNext
1704  * currentConnRate > 0 flow from iwelemNext to iwelem
1705  * With this convention, currentConnRate < 0 at the last connection for a producer
1706  * currentConnRate > 0 at the last connection for a injector
1707  */
1708 
1709  localIndex const iwelemNext = m_nextWellElemIndex[iwelem];
1710  real64 const currentConnRate = m_connRate[iwelem];
1711  localIndex iwelemUp = -1;
1712 
1713  if( iwelemNext < 0 && !m_isProducer ) // exit connection, injector
1714  {
1715  // we still need to define iwelemUp for Jacobian assembly
1716  iwelemUp = iwelem;
1717 
1718  // just copy the injection stream into compFrac
1719  for( integer ic = 0; ic < NC; ++ic )
1720  {
1721  compFracUp[ic] = m_injection[ic];
1722  for( integer jc = 0; jc < NC; ++jc )
1723  {
1724  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1725  }
1726  for( integer jc = 0; jc < NC; ++jc )
1727  {
1728  dCompFlux[ic][WJ_COFFSET::dC+jc] = 0.0;
1729  }
1730  }
1731  }
1732  else
1733  {
1734  // first set iwelemUp to the upstream cell
1735  if( ( iwelemNext < 0 && m_isProducer ) // exit connection, producer
1736  || currentConnRate < 0 ) // not an exit connection, iwelem is upstream
1737  {
1738  iwelemUp = iwelem;
1739  }
1740  else // not an exit connection, iwelemNext is upstream
1741  {
1742  iwelemUp = iwelemNext;
1743  }
1744 
1745  // copy the vars of iwelemUp into compFrac
1746  for( integer ic = 0; ic < NC; ++ic )
1747  {
1748  compFracUp[ic] = m_wellElemCompFrac[iwelemUp][ic];
1749  for( integer jc = 0; jc < NC; ++jc )
1750  {
1751  dCompFlux[ic][WJ_COFFSET::dC+jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1752  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1753  }
1754  }
1755  }
1756 
1757  // Step 2) compute upstream transport coefficient
1758 
1759  for( integer ic = 0; ic < NC; ++ic )
1760  {
1761  compFlux[ic] = compFracUp[ic] * currentConnRate;
1762  dCompFlux[ic][WJ_COFFSET::dQ] = compFracUp[ic];
1763  // none of these quantities depend on pressure
1764  dCompFlux[ic][WJ_COFFSET::dP] = 0.0;
1765  if constexpr ( IS_THERMAL )
1766  {
1767  dCompFlux[ic][WJ_COFFSET::dT] = 0.0;
1768  }
1769  for( integer jc = 0; jc < NC; ++jc )
1770  {
1771  dCompFlux[ic][WJ_COFFSET::dC+jc] = dCompFlux[ic][WJ_COFFSET::dC+jc] * currentConnRate;
1772  }
1773  }
1774 
1775  stack.offsetUp = m_wellElemDofNumber[iwelemUp];
1776  stack.iwelemUp = iwelemUp;
1777  stack.offsetCurrent = m_wellElemDofNumber[iwelem];
1778  stack.iwelemCurrent= iwelem;
1779 
1780 
1781  if( iwelemNext < 0 ) // exit connection
1782  {
1783  // for this case, we only need NC mass conservation equations
1784  computeExit ( m_dt,
1785  compFlux,
1786  stack,
1787  dCompFlux );
1788 
1789  }
1790  else // not an exit connection
1791  {
1792  compute( m_dt,
1793  compFlux,
1794  stack,
1795  dCompFlux
1796  );
1797  stack.offsetNext = m_wellElemDofNumber[iwelemNext];
1798  }
1799  stack.iwelemNext = iwelemNext;
1800  compFluxKernelOp( iwelemNext, iwelemUp, currentConnRate, dComp );
1801  }
1802 
1803 
1811  template< typename POLICY, typename KERNEL_TYPE >
1812  static void
1813  launch( localIndex const numElements,
1814  KERNEL_TYPE const & kernelComponent )
1815  {
1817  forAll< POLICY >( numElements, [=] GEOS_HOST_DEVICE ( localIndex const ie )
1818  {
1819  typename KERNEL_TYPE::StackVariables stack( 1 );
1820 
1821  kernelComponent.setup( ie, stack );
1822  kernelComponent.computeFlux( ie, stack );
1823  kernelComponent.complete( ie, stack );
1824  } );
1825  }
1826 
1827 protected:
1829  real64 const m_dt;
1832 
1837 
1840 
1841 
1846 
1851 
1854 
1856  bool const m_isProducer;
1857 
1860 
1861 
1862 };
1863 
1868 {
1869 public:
1870 
1884  template< typename POLICY >
1885  static void
1886  createAndLaunch( integer const numComps,
1887  real64 const dt,
1888  globalIndex const rankOffset,
1889  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags,
1890  string const dofKey,
1891  WellControls const & wellControls,
1892  WellElementSubRegion const & subRegion,
1893  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1894  arrayView1d< real64 > const & localRhs )
1895  {
1896  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComps, [&]( auto NC )
1897  {
1898  integer constexpr NUM_COMP = NC();
1899 
1900  using kernelType = FaceBasedAssemblyKernel< NUM_COMP, 0 >;
1901  kernelType kernel( dt, rankOffset, dofKey, wellControls, subRegion, localMatrix, localRhs, kernelFlags );
1902  kernelType::template launch< POLICY >( subRegion.size(), kernel );
1903  } );
1904  }
1905 };
1906 } // end namespace compositionalMultiphaseWellKernels
1907 
1908 } // end namespace geos
1909 
1910 #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.
arrayView4d< real64 const, constitutive::multifluid::USD_PHASE_COMP > const m_phaseCompFrac_n
Views on the phase component fraction.
GEOS_HOST_DEVICE integer elemGhostRank(localIndex const ei) const
Getter for the ghost rank of an element.
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.
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 timeAtEndOfStep, 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::SolutionCheckKernel::StackVariables createAndLaunch(integer const allowCompDensChopping, CompositionalMultiphaseFVM::ScalingType const scalingType, real64 const scalingFactor, arrayView1d< real64 const > const pressure, arrayView2d< real64 const, compflow::USD_COMP > const compDens, arrayView1d< real64 > pressureScalingFactor, arrayView1d< real64 > compDensScalingFactor, globalIndex const rankOffset, integer const numComp, string const dofKey, ElementSubRegionBase &subRegion, arrayView1d< real64 const > const localSolution)
Create a new kernel and launch.
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:1273
localIndex size() const
Get the "size" of the group, which determines the number of elements in resizable wrappers.
Definition: Group.hpp:1315
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:180
StackArray< T, 2, MAXSIZE > stackArray2d
Alias for 2D stack array.
Definition: DataTypes.hpp:204
LvArray::CRSMatrixView< T, COL_INDEX, localIndex const, LvArray::ChaiBuffer > CRSMatrixView
Alias for CRS Matrix View.
Definition: DataTypes.hpp:310
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:88
StackArray< T, 1, MAXSIZE > stackArray1d
Alias for 1D stack array.
Definition: DataTypes.hpp:188
ArraySlice< T, 2, USD > arraySlice2d
Alias for 2D array slice.
Definition: DataTypes.hpp:200
ArraySlice< T, 3, USD > arraySlice3d
Alias for 3D array slice.
Definition: DataTypes.hpp:216
ArrayView< T, 5, USD > arrayView5d
Alias for 5D array view.
Definition: DataTypes.hpp:244
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
ArraySlice< T, 1, USD > arraySlice1d
Alias for 1D array slice.
Definition: DataTypes.hpp:184
std::int32_t integer
Signed integer type.
Definition: DataTypes.hpp:82
ArrayView< T, 4, USD > arrayView4d
Alias for 4D array view.
Definition: DataTypes.hpp:228
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
ArrayView< T, 3, USD > arrayView3d
Alias for 3D array view.
Definition: DataTypes.hpp:212
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.