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 #include "physicsSolvers/fluidFlow/wells/WellPhaseVolumeRateConstraint.hpp"
44 
45 namespace geos
46 {
47 
48 
49 namespace compositionalMultiphaseWellKernels
50 {
51 
52 static constexpr real64 minDensForDivision = 1e-10;
53 
54 // tag to access well and reservoir elements in perforation rates computation
56 {
57  static constexpr integer RES = 0;
58  static constexpr integer WELL = 1;
59 };
60 
61 // tag to access the next and current well elements of a connection
62 struct ElemTag
63 {
64  static constexpr integer CURRENT = 0;
65  static constexpr integer NEXT = 1;
66 };
67 
68 // define the column offset of the derivatives
69 struct ColOffset
70 {
71  static constexpr integer DPRES = 0;
72  static constexpr integer DCOMP = 1;
73 };
74 
75 template< integer NC, integer IS_THERMAL >
77 
78 template< integer NC >
79 struct ColOffset_WellJac< NC, 0 >
80 {
81  static constexpr integer dP = 0;
82  static constexpr integer dC = 1;
83  static constexpr integer dQ = dC + NC;
84  static integer constexpr nDer = dQ + 1;
85 
86 };
87 
88 template< integer NC >
89 struct ColOffset_WellJac< NC, 1 >
90 {
91  static constexpr integer dP = 0;
92  static constexpr integer dC = 1;
93  static constexpr integer dQ = dC + NC;
94  static constexpr integer dT = dQ+1;
96  static integer constexpr nDer = dT + 1;
97 };
98 
99 // define the row offset of the residual equations
100 struct RowOffset
101 {
102  static constexpr integer CONTROL = 0;
103  static constexpr integer MASSBAL = 1;
104 };
105 
106 template< integer NC, integer IS_THERMAL >
108 
109 template< integer NC >
110 struct RowOffset_WellJac< NC, 0 >
111 {
112  static constexpr integer CONTROL = 0;
113  static constexpr integer MASSBAL = 1;
114  static constexpr integer VOLBAL = MASSBAL + NC;
115  static constexpr integer nEqn = VOLBAL+1;
116 };
117 
118 template< integer NC >
119 struct RowOffset_WellJac< NC, 1 >
120 {
121  static constexpr integer CONTROL = 0;
122  static constexpr integer MASSBAL = 1;
123  static constexpr integer VOLBAL = MASSBAL + NC;
124  static constexpr integer ENERGYBAL = VOLBAL+1;
125  static constexpr integer nEqn = ENERGYBAL+1;
126 
127 };
128 /******************************** ControlEquationHelper ********************************/
130 {
133 
135  inline
136  static
137  void
138  selectLimitingConstraint( bool const isProducer,
139  WellControls::Control const & inputControl,
140  WellControls::Control const & currentControl,
141  integer const phasePhaseIndex,
142  real64 const & targetBHP,
143  real64 const & targetPhaseRate,
144  real64 const & targetTotalRate,
145  real64 const & targetMassRate,
146  real64 const & currentBHP,
147  arrayView1d< real64 const > const & currentPhaseVolRate,
148  real64 const & currentTotalVolRate,
149  real64 const & currentMassRate,
150  WellControls::Control & newControl );
151 
152  template< integer NC, integer IS_THERMAL >
154  inline
155  static void
156  compute( globalIndex const rankOffset,
157  WellControls::Control const currentControl,
158  integer const targetPhaseIndex,
159  real64 const & targetBHP,
160  real64 const & targetPhaseRate,
161  real64 const & targetTotalRate,
162  real64 const & targetMassRate,
163  real64 const & currentBHP,
164  real64 const & targetValue,
165  arrayView1d< real64 const > const & dCurrentBHP,
166  arrayView1d< real64 const > const & currentPhaseVolRate,
167  arrayView2d< real64 const > const & dCurrentPhaseVolRate,
168  real64 const & currentTotalVolRate,
169  arrayView1d< real64 const > const & dCurrentTotalVolRate,
170  real64 const & massDensity,
171  globalIndex const dofNumber,
172  CRSMatrixView< real64, globalIndex const > const & localMatrix,
173  arrayView1d< real64 > const & localRhs );
174 
175 };
176 
177 /******************************** PressureRelationKernel ********************************/
178 
180 {
181  using Deriv = constitutive::multifluid::DerivativeOffset;
185 
186  template< integer NC, integer IS_THERMAL >
188  inline
189  static void
190  compute( real64 const & gravCoef,
191  real64 const & gravCoefNext,
192  real64 const & pres,
193  real64 const & presNext,
194  real64 const & totalMassDens,
195  real64 const & totalMassDensNext,
198  real64 & localPresRel,
199  real64 ( &localPresRelJacobian )[2*(NC+1+IS_THERMAL)] );
200 
201  template< integer NC, integer IS_THERMAL >
202  static void
203  launch( localIndex const size,
204  globalIndex const rankOffset,
205  arrayView1d< integer const > const elemStatus,
206  arrayView1d< globalIndex const > const & wellElemDofNumber,
207  arrayView1d< real64 const > const & wellElemGravCoef,
208  arrayView1d< localIndex const > const & nextWellElemIndex,
209  arrayView1d< real64 const > const & wellElemPressure,
210  arrayView1d< real64 const > const & wellElemTotalMassDens,
211  arrayView2d< real64 const, compflow::USD_FLUID_DC > const & dWellElemTotalMassDens,
212  bool & controlHasSwitched,
213  CRSMatrixView< real64, globalIndex const > const & localMatrix,
214  arrayView1d< real64 > const & localRhs );
215 
216 };
217 
218 /******************************** VolumeBalanceKernel ********************************/
219 
221 {
222 
225 
226  template< integer NC >
228  inline
229  static void
230  compute( integer const numPhases,
231  real64 const & volume,
234  real64 & localVolBalance,
235  real64 ( &localVolBalanceJacobian )[NC+1] );
236 
237  template< integer NC >
238  static void
239  launch( localIndex const size,
240  integer const numPhases,
241  globalIndex const rankOffset,
242  arrayView1d< globalIndex const > const & wellElemDofNumber,
243  arrayView1d< integer const > const & wellElemGhostRank,
244  arrayView2d< real64 const, compflow::USD_PHASE > const & wellElemPhaseVolFrac,
245  arrayView3d< real64 const, compflow::USD_PHASE_DC > const & dWellElemPhaseVolFrac,
246  arrayView1d< real64 const > const & wellElemVolume,
247  CRSMatrixView< real64, globalIndex const > const & localMatrix,
248  arrayView1d< real64 > const & localRhs );
249 
250 };
251 
252 /******************************** PresTempCompFracInitializationKernel ********************************/
253 
255 {
256 
257  using CompFlowAccessors =
258  StencilAccessors< fields::flow::pressure,
259  fields::flow::temperature,
260  fields::flow::globalCompDensity,
261  fields::flow::phaseVolumeFraction >;
262 
263  using MultiFluidAccessors =
264  StencilMaterialAccessors< constitutive::MultiFluidBase,
265  fields::multifluid::phaseMassDensity >;
266 
267 
275  template< typename VIEWTYPE >
277 
278  static void
279  launch( localIndex const perforationSize,
280  localIndex const subRegionSize,
281  integer const numComponents,
282  integer const numPhases,
283  WellControls const & wellControls,
284  real64 const & currentTime,
290  arrayView1d< localIndex const > const & resElementRegion,
291  arrayView1d< localIndex const > const & resElementSubRegion,
292  arrayView1d< localIndex const > const & resElementIndex,
293  arrayView1d< real64 const > const & perfGravCoef,
294  arrayView1d< integer const > const & perfState,
295  arrayView1d< real64 const > const & wellElemGravCoef,
296  arrayView1d< real64 > const & wellElemPres,
297  arrayView1d< real64 > const & wellElemTemp,
298  arrayView2d< real64, compflow::USD_COMP > const & wellElemCompFrac );
299 
300 };
301 
302 /******************************** CompDensInitializationKernel ********************************/
303 
305 {
306 
307  static void
308  launch( localIndex const subRegionSize,
309  integer const numComponents,
310  arrayView2d< real64 const, compflow::USD_COMP > const & wellElemCompFrac,
312  arrayView2d< real64, compflow::USD_COMP > const & wellElemCompDens );
313 
314 };
315 
316 /******************************** RateInitializationKernel ********************************/
317 
319 {
320 
321  static void
322  launch( localIndex const subRegionSize,
323  WellControls const & wellControls,
324  real64 const & currentTime,
327  arrayView1d< real64 > const & connRate );
328 
329 };
330 
331 
332 /******************************** TotalMassDensityKernel ****************************/
333 
340 template< integer NUM_COMP, integer NUM_PHASE >
342 {
343 public:
344 
346  using Base::numComp;
347 
349  static constexpr integer numPhase = NUM_PHASE;
350 
357  constitutive::MultiFluidBase const & fluid )
358  : Base(),
359  m_phaseVolFrac( subRegion.getField< fields::well::phaseVolumeFraction >() ),
360  m_dPhaseVolFrac( subRegion.getField< fields::well::dPhaseVolumeFraction >() ),
361  m_dCompFrac_dCompDens( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
362  m_phaseMassDens( fluid.phaseMassDensity() ),
363  m_dPhaseMassDens( fluid.dPhaseMassDensity() ),
364  m_totalMassDens( subRegion.getField< fields::well::totalMassDensity >() ),
365  m_dTotalMassDens( subRegion.getField< fields::well::dTotalMassDensity >() )
366  {}
367 
374  template< typename FUNC = NoOpFunc >
376  inline
377  void compute( localIndex const ei,
378  FUNC && totalMassDensityKernelOp = NoOpFunc{} ) const
379  {
380  using Deriv = constitutive::multifluid::DerivativeOffset;
381 
382  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
383  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
384  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
385  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseMassDens = m_phaseMassDens[ei][0];
386  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseMassDens = m_dPhaseMassDens[ei][0];
387  real64 & totalMassDens = m_totalMassDens[ei];
388  arraySlice1d< real64, compflow::USD_FLUID_DC - 1 > dTotalMassDens = m_dTotalMassDens[ei];
389 
390  real64 dMassDens_dC[numComp]{};
391 
392  totalMassDens = 0.0;
393 
394  dTotalMassDens[Deriv::dP]=0.0;
395  for( integer ic = 0; ic < numComp; ++ic )
396  {
397  dTotalMassDens[Deriv::dC+ic]=0.0;
398  }
399 
400  for( integer ip = 0; ip < numPhase; ++ip )
401  {
402  totalMassDens += phaseVolFrac[ip] * phaseMassDens[ip];
403  dTotalMassDens[Deriv::dP] += dPhaseVolFrac[ip][Deriv::dP] * phaseMassDens[ip] + phaseVolFrac[ip] * dPhaseMassDens[ip][Deriv::dP];
404 
405  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseMassDens[ip], dMassDens_dC, Deriv::dC );
406  for( integer ic = 0; ic < numComp; ++ic )
407  {
408  dTotalMassDens[Deriv::dC+ic] += dPhaseVolFrac[ip][Deriv::dC+ic] * phaseMassDens[ip]
409  + phaseVolFrac[ip] * dMassDens_dC[ic];
410  }
411 
412  totalMassDensityKernelOp( ip ); //, phaseVolFrac, dTotalMassDens_dPres, dTotalMassDens_dCompDens );
413  }
414 
415  }
416 
417 protected:
418 
419  // inputs
420 
425 
429 
430  // outputs
431 
435 
436 
437 };
438 
443 {
444 public:
445 
454  template< typename POLICY >
455  static void
456  createAndLaunch( integer const numComp,
457  integer const numPhase,
458  ObjectManagerBase & subRegion,
459  constitutive::MultiFluidBase const & fluid )
460  {
461  if( numPhase == 2 )
462  {
463  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
464  {
465  integer constexpr NUM_COMP = NC();
466  TotalMassDensityKernel< NUM_COMP, 2 > kernel( subRegion, fluid );
467  TotalMassDensityKernel< NUM_COMP, 2 >::template launch< POLICY >( subRegion.size(), kernel );
468  } );
469  }
470  else if( numPhase == 3 )
471  {
472  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
473  {
474  integer constexpr NUM_COMP = NC();
475  TotalMassDensityKernel< NUM_COMP, 3 > kernel( subRegion, fluid );
476  TotalMassDensityKernel< NUM_COMP, 3 >::template launch< POLICY >( subRegion.size(), kernel );
477  } );
478  }
479  }
480 };
481 
482 
483 /******************************** ResidualNormKernel ********************************/
484 
489 {
490 public:
491 
493  using Base::m_minNormalizer;
494  using Base::m_rankOffset;
495  using Base::m_localResidual;
496  using Base::m_dofNumber;
497 
498  ResidualNormKernel( globalIndex const rankOffset,
499  arrayView1d< real64 const > const & localResidual,
500  arrayView1d< globalIndex const > const & dofNumber,
502  integer const numComp,
503  integer const numDof,
504  WellElementSubRegion const & subRegion,
505  constitutive::MultiFluidBase const & fluid,
506  WellControls const & wellControls,
507  real64 const time,
508  real64 const dt,
509  real64 const minNormalizer )
510  : Base( rankOffset,
511  localResidual,
512  dofNumber,
513  ghostRank,
514  minNormalizer ),
515  m_numComp( numComp ),
516  m_numDof( numDof ),
517  m_dt( dt ),
518  m_isLocallyOwned( subRegion.isLocallyOwned() ),
520  m_isProducer( wellControls.isProducer() ),
521  m_currentControl( wellControls.getControl() ),
522  m_volume( subRegion.getElementVolume() ),
523  m_phaseDens_n( fluid.phaseDensity_n() ),
524  m_totalDens_n( fluid.totalDensity_n() )
525  {
526  if( m_isProducer )
527  {
528  m_targetBHP = wellControls.getMinBHPConstraint()->getConstraintValue( time );
529  // Note this assumes that there is only one rate constraint
530  // This is a normalizer for the balance equations. The normalizaer should be the current rate not the constraint value!!
531  // This is one of the reasons for restricting constraint type for a production well
532  // another pr will remove fix this (so the cause for difference results is isolated to one change)
534  m_constraintValue = wellControls.getProdRateConstraints()[0]->getConstraintValue( time );
535 
536  }
537  else
538  {
539  m_targetBHP = wellControls.getMaxBHPConstraint()->getConstraintValue( time );
540 
541  // Note this assumes that there is only one rate constraint
542  // This is a normalizer for the balance equations. The normalizaer should be the current rate not the constraint value!!
543  // This is one of the reasons for restricting constraint type for a production well
544  // another pr will remove fix this (so the cause for difference results is isolated to one change)
545  m_targetPhaseIndex = -1;
546  m_constraintValue = wellControls.getInjRateConstraints()[0]->getConstraintValue( time );
547 
548  }
549 
550 
551  }
552 
554  virtual void computeLinf( localIndex const iwelem,
555  LinfStackVariables & stack ) const override
556  {
558 
559  real64 normalizer = 0.0;
560  for( integer idof = 0; idof < m_numDof; ++idof )
561  {
562 
563  // Step 1: compute a normalizer for the control or pressure equation
564 
565  // for the control equation, we distinguish two cases
566  if( idof == ROFFSET::CONTROL )
567  {
568 
569  // for the top well element, normalize using the current control
570  if( m_isLocallyOwned && iwelem == m_iwelemControl )
571  {
573  {
574  // the residual entry is in pressure units
575  normalizer = m_targetBHP;
576  }
578  {
579  // the residual entry is in volume / time units
580  normalizer = LvArray::math::max( LvArray::math::abs( m_constraintValue ), m_minNormalizer );
581  }
583  {
584  // the residual entry is in volume / time units
585  normalizer = LvArray::math::max( LvArray::math::abs( m_constraintValue ), m_minNormalizer );
586  }
588  {
589  // the residual entry is in volume / time units
590  normalizer = LvArray::math::max( LvArray::math::abs( m_constraintValue ), m_minNormalizer );
591  }
592  }
593  // for the pressure difference equation, always normalize by the BHP
594  else
595  {
596  normalizer = m_targetBHP;
597  }
598  }
599  // Step 2: compute a normalizer for the mass balance equations
600  else if( idof >= ROFFSET::MASSBAL && idof < ROFFSET::MASSBAL + m_numComp )
601  {
602  if( m_isProducer ) // only PHASEVOLRATE is supported for now
603  {
604  // the residual is in mass units
605  normalizer = m_dt * LvArray::math::abs( m_constraintValue ) * m_phaseDens_n[iwelem][0][m_targetPhaseIndex];
606  }
607  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
608  {
610  {
611  normalizer = m_dt * LvArray::math::abs( m_constraintValue );
612  }
613  else
614  {
615  // the residual is in mass units
616  normalizer = m_dt * LvArray::math::abs( m_constraintValue ) * m_totalDens_n[iwelem][0];
617  }
618 
619  }
620 
621  // to make sure that everything still works well if the rate is zero, we add this check
622  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] * m_totalDens_n[iwelem][0] );
623  }
624  // Step 3: compute a normalizer for the volume balance equations
625  else if( idof == ROFFSET::MASSBAL + m_numComp )
626  {
627  if( m_isProducer ) // only PHASEVOLRATE is supported for now
628  {
629  // the residual is in volume units
630  normalizer = m_dt * LvArray::math::abs( m_constraintValue );
631  }
632  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
633  {
635  {
636  normalizer = m_dt * LvArray::math::abs( m_constraintValue/ m_totalDens_n[iwelem][0] );
637  }
638  else
639  {
640  normalizer = m_dt * LvArray::math::abs( m_constraintValue );
641  }
642 
643  }
644 
645  }
646 
647  // to make sure that everything still works well if the rate is zero, we add this check
648  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] );
649 
650  // Step 4: compute the contribution to the residual
651  real64 const val = LvArray::math::abs( m_localResidual[stack.localRow + idof] ) / normalizer;
652  if( val > stack.localValue[0] )
653  {
654  stack.localValue[0] = val;
655  }
656  }
657  }
658 
660  virtual void computeL2( localIndex const iwelem,
661  L2StackVariables & stack ) const override
662  {
663  GEOS_UNUSED_VAR( iwelem, stack );
664  GEOS_ERROR( "The L2 norm is not implemented for CompositionalMultiphaseWell" );
665  }
666 
667 
668 protected:
669 
672 
675 
678 
680  real64 const m_dt;
681 
683  bool const m_isLocallyOwned;
684 
687 
689  bool const m_isProducer;
690 
693  real64 m_constraintValue;
694  real64 m_targetBHP;
695 
696 
699 
703 
704 };
705 
710 {
711 public:
712 
728  template< typename POLICY >
729  static void
730  createAndLaunch( integer const numComp,
731  integer const numDof,
732  globalIndex const rankOffset,
733  string const & dofKey,
734  arrayView1d< real64 const > const & localResidual,
735  WellElementSubRegion const & subRegion,
736  constitutive::MultiFluidBase const & fluid,
737  WellControls const & wellControls,
738  real64 const time,
739  real64 const dt,
740  real64 const minNormalizer,
741  real64 (& residualNorm)[1] )
742  {
743  arrayView1d< globalIndex const > const dofNumber = subRegion.getReference< array1d< globalIndex > >( dofKey );
744  arrayView1d< integer const > const ghostRank = subRegion.ghostRank();
745 
746  ResidualNormKernel kernel( rankOffset, localResidual, dofNumber, ghostRank,
747  numComp, numDof, subRegion, fluid, wellControls, time, dt, minNormalizer );
748  ResidualNormKernel::launchLinf< POLICY >( subRegion.size(), kernel, residualNorm );
749  }
750 
751 };
752 
753 /******************************** SolutionScalingKernel ********************************/
754 
759 {
760 public:
761 
775  template< typename POLICY >
777  createAndLaunch( real64 const maxRelativePresChange,
778  real64 const maxAbsolutePresChange,
779  real64 const maxCompFracChange,
780  real64 const maxRelativeCompDensChange,
781  globalIndex const rankOffset,
782  integer const numComp,
783  string const dofKey,
784  ElementSubRegionBase & subRegion,
785  arrayView1d< real64 const > const localSolution )
786  {
787  arrayView1d< real64 const > const pressure =
788  subRegion.getField< fields::well::pressure >();
790  subRegion.getField< fields::well::globalCompDensity >();
791  arrayView1d< real64 > pressureScalingFactor =
792  subRegion.getField< fields::well::pressureScalingFactor >();
793  arrayView1d< real64 > compDensScalingFactor =
794  subRegion.getField< fields::well::globalCompDensityScalingFactor >();
795  isothermalCompositionalMultiphaseBaseKernels::
796  SolutionScalingKernel kernel( maxRelativePresChange, maxAbsolutePresChange, maxCompFracChange, maxRelativeCompDensChange, rankOffset,
797  numComp, dofKey, subRegion, localSolution, pressure, compDens, pressureScalingFactor, compDensScalingFactor );
798  return isothermalCompositionalMultiphaseBaseKernels::
799  SolutionScalingKernel::
800  launch< POLICY >( subRegion.size(), kernel );
801  }
802 
803 };
804 
805 /******************************** ElementBasedAssemblyKernel ********************************/
806 
813 template< integer NUM_COMP, integer IS_THERMAL >
815 {
816 public:
819 
820  // Well jacobian column and row indicies
821  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NUM_COMP, IS_THERMAL >;
825  static constexpr integer numComp = NUM_COMP;
826 
828  static constexpr integer numDof = NUM_COMP + 1 + IS_THERMAL;
829 
831  static constexpr integer numEqn = NUM_COMP + 1 + IS_THERMAL;
832 
833 
846  integer const isProducer,
847  globalIndex const rankOffset,
848  string const dofKey,
849  WellElementSubRegion const & subRegion,
850  constitutive::MultiFluidBase const & fluid,
851  CRSMatrixView< real64, globalIndex const > const & localMatrix,
852  arrayView1d< real64 > const & localRhs,
853  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const kernelFlags )
854  : m_numPhases( numPhases ),
855  m_isProducer( isProducer ),
856  m_rankOffset( rankOffset ),
857  m_iwelemControl( subRegion.getTopWellElementIndex() ),
858  m_dofNumber( subRegion.getReference< array1d< globalIndex > >( dofKey ) ),
859  m_elemGhostRank( subRegion.ghostRank() ),
860  m_elemStatus( subRegion.getLocalWellElementStatus() ),
861  m_volume( subRegion.getElementVolume() ),
862  m_dCompFrac_dCompDens( subRegion.getField< fields::flow::dGlobalCompFraction_dGlobalCompDensity >() ),
863  m_phaseVolFrac_n( subRegion.getField< fields::flow::phaseVolumeFraction_n >() ),
864  m_phaseVolFrac( subRegion.getField< fields::flow::phaseVolumeFraction >() ),
865  m_dPhaseVolFrac( subRegion.getField< fields::flow::dPhaseVolumeFraction >() ),
866  m_phaseDens_n( fluid.phaseDensity_n() ),
867  m_phaseDens( fluid.phaseDensity() ),
868  m_dPhaseDens( fluid.dPhaseDensity() ),
869  m_phaseCompFrac_n( fluid.phaseCompFraction_n() ),
870  m_phaseCompFrac( fluid.phaseCompFraction() ),
871  m_dPhaseCompFrac( fluid.dPhaseCompFraction() ),
872  m_compDens( subRegion.getField< fields::flow::globalCompDensity >() ),
873  m_compDens_n( subRegion.getField< fields::flow::globalCompDensity_n >() ),
874  m_localMatrix( localMatrix ),
875  m_localRhs( localRhs ),
876  m_kernelFlags( kernelFlags )
877  {}
878 
884  {
885 public:
886 
887  // volume information (used by both accumulation and volume balance)
888  real64 volume = 0.0;
889 
890  // Residual information
891 
894 
896  globalIndex dofIndices[numDof]{}; // NC compdens + P + thermal
897  globalIndex eqnRowIndices[numDof]{};
898  globalIndex dofColIndices[numDof]{};
899 
902 
905 
906  };
907 
914  bool skipElement( localIndex const ei ) const
915  {
916  return ( m_elemGhostRank( ei ) >= 0 || m_elemStatus[ei]==WellElementSubRegion::WellElemStatus::CLOSED );
917  }
918 
925  void setup( localIndex const ei,
926  StackVariables & stack ) const
927  {
928  // initialize the volume
929  stack.volume = m_volume[ei];
930 
931  // Note row/col indices needed to be consistent with layout of stack.localJacobian
932  // Setup row equation indices for this element ( mass + vol + thermal if valid)
933 
934  // 1) Mass Balance
935  for( integer ic = 0; ic < numComp; ++ic )
936  {
937  stack.eqnRowIndices[ic] = m_dofNumber[ei] + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
938  }
939 // 2) Volume Balance
940  stack.eqnRowIndices[numComp] = m_dofNumber[ei] + WJ_ROFFSET::VOLBAL - m_rankOffset;
941  // 3) Energy Balance
942  if constexpr ( IS_THERMAL )
943  {
944  stack.eqnRowIndices[numComp+1] = m_dofNumber[ei] + WJ_ROFFSET::ENERGYBAL - m_rankOffset;
945  }
946  // Setup equation column indices for this element ( P + COMPDENS + THERMAL if valid)
947  stack.dofColIndices[0] = m_dofNumber[ei] + WJ_COFFSET::dP;
948  for( integer ic = 0; ic < numComp; ++ic )
949  {
950  stack.dofColIndices[ic+1] = m_dofNumber[ei] + WJ_COFFSET::dC+ic;
951  }
952  if constexpr ( IS_THERMAL )
953  {
954  stack.dofColIndices[numComp+1] = m_dofNumber[ei] + WJ_COFFSET::dT;
955  }
956  for( integer jc = 0; jc < numEqn; ++jc )
957  {
958  stack.localResidual[jc] = 0.0;
959  for( integer ic = 0; ic < numDof; ++ic )
960  {
961  stack.localJacobian[jc][ic] = 0.0;
962  }
963  }
964  }
972  template< typename FUNC = NoOpFunc >
975  StackVariables & stack,
976  FUNC && phaseAmountKernelOp = NoOpFunc{} ) const
977  {
978 
979  using Deriv = constitutive::multifluid::DerivativeOffset;
980 
981  // construct the slices for variables accessed multiple times
982  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
983 
984  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac_n = m_phaseVolFrac_n[ei];
985  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
986  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
987 
988  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens_n = m_phaseDens_n[ei][0];
989  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens = m_phaseDens[ei][0];
990  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseDens = m_dPhaseDens[ei][0];
991 
992  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac_n = m_phaseCompFrac_n[ei][0];
993  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac = m_phaseCompFrac[ei][0];
994  arraySlice3d< real64 const, constitutive::multifluid::USD_PHASE_COMP_DC - 2 > dPhaseCompFrac = m_dPhaseCompFrac[ei][0];
995 
996  // temporary work arrays
997  real64 dPhaseAmount[FLUID_PROP_COFFSET::nDer]{};
998  real64 dPhaseAmount_dC[numComp]{};
999  real64 dPhaseCompFrac_dC[numComp]{};
1000 
1001  // sum contributions to component accumulation from each phase
1002  for( integer ip = 0; ip < m_numPhases; ++ip )
1003  {
1004  real64 const phaseAmount = stack.volume * phaseVolFrac[ip] * phaseDens[ip];
1005  real64 const phaseAmount_n = stack.volume * phaseVolFrac_n[ip] * phaseDens_n[ip];
1006 
1007  dPhaseAmount[FLUID_PROP_COFFSET::dP]=stack.volume * ( dPhaseVolFrac[ip][Deriv::dP] * phaseDens[ip]
1008  + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dP] );
1009 
1010  // assemble density dependence
1011  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], dPhaseAmount_dC, Deriv::dC );
1012  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], &dPhaseAmount[FLUID_PROP_COFFSET::dC], Deriv::dC );
1013  for( integer jc = 0; jc < numComp; ++jc )
1014  {
1015  dPhaseAmount_dC[jc] = dPhaseAmount_dC[jc] * phaseVolFrac[ip]
1016  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1017  dPhaseAmount_dC[jc] *= stack.volume;
1018  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] = dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] * phaseVolFrac[ip]
1019  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1020  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] *= stack.volume;
1021  }
1022  // ic - index of component whose conservation equation is assembled
1023  // (i.e. row number in local matrix)
1024  for( integer ic = 0; ic < numComp; ++ic )
1025  {
1026  real64 const phaseCompAmount = phaseAmount * phaseCompFrac[ip][ic];
1027  real64 const phaseCompAmount_n = phaseAmount_n * phaseCompFrac_n[ip][ic];
1028 
1029  real64 const dPhaseCompAmount_dP = dPhaseAmount[FLUID_PROP_COFFSET::dP] * phaseCompFrac[ip][ic]
1030  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dP];
1031 
1032  stack.localResidual[ic] += phaseCompAmount - phaseCompAmount_n;
1033  stack.localJacobian[ic][0] += dPhaseCompAmount_dP;
1034 
1035  // jc - index of component w.r.t. whose compositional var the derivative is being taken
1036  // (i.e. col number in local matrix)
1037 
1038  // assemble phase composition dependence
1039  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseCompFrac[ip][ic], dPhaseCompFrac_dC, Deriv::dC );
1040  for( integer jc = 0; jc < numComp; ++jc )
1041  {
1042  real64 const dPhaseCompAmount_dC = dPhaseCompFrac_dC[jc] * phaseAmount
1043  + phaseCompFrac[ip][ic] * dPhaseAmount[FLUID_PROP_COFFSET::dC+jc];
1044 
1045  stack.localJacobian[ic][jc + 1] += dPhaseCompAmount_dC;
1046  }
1047  }
1048  if constexpr ( IS_THERMAL )
1049  {
1050  dPhaseAmount[FLUID_PROP_COFFSET::dT] = stack.volume * (dPhaseVolFrac[ip][Deriv::dT] * phaseDens[ip] + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dT] );
1051  for( integer ic = 0; ic < numComp; ++ic )
1052  {
1053  // assemble the derivatives of the component mass balance equations with respect to temperature
1054  stack.localJacobian[ic][numComp+1] += dPhaseAmount[FLUID_PROP_COFFSET::dT] * phaseCompFrac[ip][ic]
1055  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dT];
1056  }
1057  }
1058  // call the lambda in the phase loop to allow the reuse of the phase amounts and their derivatives
1059  // possible use: assemble accumulation term of the energy equation for this phase
1060  phaseAmountKernelOp( ip, phaseAmount, phaseAmount_n, dPhaseAmount );
1061 
1062  }
1063 
1064  // check zero diagonal (works only in debug)
1065  /*
1066  for( integer ic = 0; ic < numComp; ++ic )
1067  {
1068  GEOS_ASSERT_MSG ( LvArray::math::abs( stack.localJacobian[ic][ic] ) > minDensForDivision,
1069  GEOS_FMT( "Zero diagonal in Jacobian: equation {}, value = {}", ic, stack.localJacobian[ic][ic] ) );
1070  }
1071  */
1072  }
1073 
1074 
1085  StackVariables & stack ) const
1086  {
1087  using Deriv = constitutive::multifluid::DerivativeOffset;
1088 
1089  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
1090  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
1091 
1092  real64 oneMinusPhaseVolFracSum = 1.0;
1093 
1094  // sum contributions to component accumulation from each phase
1095 // Note localJacobian stores equation balances in order of component/vol/enerqy
1096  // These are mapped to solver orderings with indicies setup in stack variables
1097  for( integer ip = 0; ip < m_numPhases; ++ip )
1098  {
1099  oneMinusPhaseVolFracSum -= phaseVolFrac[ip];
1100  stack.localJacobian[numComp][0] -= dPhaseVolFrac[ip][Deriv::dP];
1101 
1102  for( integer jc = 0; jc < numComp; ++jc )
1103  {
1104  stack.localJacobian[numComp][jc+1] -= dPhaseVolFrac[ip][Deriv::dC+jc];
1105  }
1106 
1107  if constexpr ( IS_THERMAL)
1108  {
1109  stack.localJacobian[numComp][numComp+1] -= dPhaseVolFrac[ip][Deriv::dT];
1110  }
1111 
1112  }
1113  // scale saturation-based volume balance by pore volume (for better scaling w.r.t. other equations)
1114  stack.localResidual[numComp] = stack.volume * oneMinusPhaseVolFracSum;
1115  for( integer idof = 0; idof < numComp+1+IS_THERMAL; ++idof )
1116  {
1117  stack.localJacobian[numComp][idof] *= stack.volume;
1118  }
1119  }
1120 
1127  void complete( localIndex const ei, //GEOS_UNUSED_PARAM( ei ),
1128  StackVariables & stack ) const
1129  {
1130  using namespace compositionalMultiphaseUtilities;
1131 
1132  integer const numRows = numComp+1+ IS_THERMAL;
1133 
1134  if constexpr ( IS_THERMAL)
1135  {
1136  if( ei == m_iwelemControl && !m_isProducer )
1137  {
1138  // For top segment energy balance eqn replaced with T(n+1) - T = 0
1139  // No other energy balance derivatives
1140  // Assumption is global index == 0 is top segment with fixed temp BC
1141 
1142  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1143  {
1144  stack.localJacobian[numRows-1][i] = 0.0;
1145  }
1146  // constant Temperature
1147  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1148  stack.localJacobian[i][numRows-1] = 0.0;
1149  stack.localJacobian[numRows-1][numRows-1] = 1.0;
1150 
1151  stack.localResidual[numRows-1]=0.0;
1152  }
1153  }
1154 
1155  if( m_kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) )
1156  {
1157  // apply equation/variable change transformation to the component mass balance equations
1158  real64 work[numComp + 1 + IS_THERMAL]{};
1159  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numComp, numComp+1+ IS_THERMAL, stack.localJacobian, work );
1160  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numComp, stack.localResidual );
1161  }
1162 
1163  // add contribution to residual and jacobian into:
1164  // - the component mass balance equations (i = 0 to i = numComp-1)
1165  // - the volume balance equations (i = numComp)
1166  // note that numDof includes derivatives wrt temperature if this class is derived in ThermalKernels
1167 
1168  for( integer i = 0; i < numRows; ++i )
1169  {
1170  m_localRhs[stack.eqnRowIndices[i]] += stack.localResidual[i];
1171  m_localMatrix.addToRow< serialAtomic >( stack.eqnRowIndices[i],
1172  stack.dofColIndices,
1173  stack.localJacobian[i],
1174  numComp+1+ IS_THERMAL );
1175  }
1176 
1177  }
1178 
1186  template< typename POLICY, typename KERNEL_TYPE >
1187  static void
1188  launch( localIndex const numElems,
1189  KERNEL_TYPE const & kernelComponent )
1190  {
1192  forAll< POLICY >( numElems, [=] GEOS_HOST_DEVICE ( localIndex const ei )
1193  {
1194  if( kernelComponent.skipElement( ei ) )
1195  {
1196  return;
1197  }
1198 
1199  typename KERNEL_TYPE::StackVariables stack;
1200 
1201  kernelComponent.setup( ei, stack );
1202  kernelComponent.computeAccumulation( ei, stack );
1203  kernelComponent.computeVolumeBalance( ei, stack );
1204  kernelComponent.complete( ei, stack );
1205  } );
1206  }
1207 
1208 protected:
1209 
1212 
1215 
1218 
1221 
1224 
1227 
1230 
1233 
1236  arrayView2d< real64 const > const m_porosity;
1237  arrayView2d< real64 const > const m_dPoro_dPres;
1238 
1241 
1246 
1251 
1256 
1257  // Views on component densities
1260 
1265 
1266  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const m_kernelFlags;
1267 };
1268 
1269 
1274 {
1275 public:
1288  template< typename POLICY >
1289  static void
1290  createAndLaunch( localIndex const numComps,
1291  localIndex const numPhases,
1292  integer const isProducer,
1293  globalIndex const rankOffset,
1294  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags,
1295  string const dofKey,
1296  WellElementSubRegion const & subRegion,
1297  constitutive::MultiFluidBase const & fluid,
1298  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1299  arrayView1d< real64 > const & localRhs )
1300  {
1301  geos::internal::kernelLaunchSelectorCompThermSwitch( numComps, 0, [&]( auto NC, auto IS_THERMAL )
1302  {
1303  localIndex constexpr NUM_COMP = NC();
1304 
1305  integer constexpr istherm = IS_THERMAL();
1306 
1308  kernel( numPhases, isProducer, rankOffset, dofKey, subRegion, fluid, localMatrix, localRhs, kernelFlags );
1310  launch< POLICY, ElementBasedAssemblyKernel< NUM_COMP, istherm > >( subRegion.size(), kernel );
1311  } );
1312  }
1313 };
1322 template< integer NC, integer IS_THERMAL >
1324 {
1325 public:
1326 
1330 
1331  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1334 
1335  using CP_Deriv = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1337  static constexpr integer numComp = NC;
1338 
1340  static constexpr integer numDof = WJ_COFFSET::nDer;
1341 
1343  static constexpr integer numEqn = NC;
1344 
1345  static constexpr integer maxNumElems = 2;
1346  static constexpr integer maxStencilSize = 2;
1362  globalIndex const rankOffset,
1363  string const wellDofKey,
1364  WellControls const & wellControls,
1365  WellElementSubRegion const & subRegion,
1366  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1367  arrayView1d< real64 > const & localRhs,
1368  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags )
1369  :
1370  m_dt( dt ),
1371  m_rankOffset( rankOffset ),
1372  m_wellElemDofNumber ( subRegion.getReference< array1d< globalIndex > >( wellDofKey ) ),
1373  m_nextWellElemIndex ( subRegion.getReference< array1d< localIndex > >( WellElementSubRegion::viewKeyStruct::nextWellElementIndexString()) ),
1374  m_elemStatus( subRegion.getLocalWellElementStatus() ),
1375  m_connRate ( subRegion.getField< fields::well::mixtureConnectionRate >() ),
1376  m_wellElemCompFrac ( subRegion.getField< fields::well::globalCompFraction >() ),
1377  m_dWellElemCompFrac_dCompDens ( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
1378  m_localMatrix( localMatrix ),
1379  m_localRhs ( localRhs ),
1380  m_useTotalMassEquation ( kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) ),
1381  m_isProducer ( wellControls.isProducer() ),
1382  m_injection ( wellControls.getInjectionStream() )
1383  {}
1384 
1386  {
1387 public:
1388 
1396  : stencilSize( size ),
1397  numConnectedElems( 2 ),
1398  dofColIndices( size * numDof )
1399  {}
1400 
1401  // Stencil information
1402  localIndex const stencilSize;
1405 
1406 
1407  // edge indexes
1408  localIndex iwelemUp;
1409  localIndex iwelemNext;
1410  localIndex iwelemCurrent;
1411  globalIndex offsetUp;
1412  globalIndex offsetCurrent;
1413  globalIndex offsetNext;
1414  // Local degrees of freedom and local residual/jacobian
1415 
1418 
1425  };
1426 
1427 
1429  bool skipElement( localIndex const ei ) const
1430  {
1431  return ( m_elemStatus[ei]==WellElementSubRegion::WellElemStatus::CLOSED );
1432  }
1433 
1440  inline
1441  void setup( localIndex const iconn,
1442  StackVariables & stack ) const
1443  {
1444  stack.numConnectedElems=2;
1445  if( m_nextWellElemIndex[iconn] <0 )
1446  {
1447  stack.numConnectedElems = 1;
1448  }
1449 
1450  stack.localFlux.resize( stack.numConnectedElems*numEqn );
1451  stack.localFluxJacobian.resize( stack.numConnectedElems * numEqn, stack.stencilSize * numDof );
1452  stack.localFluxJacobian_dQ.resize( stack.numConnectedElems * numEqn, 1 );
1453 
1454  }
1455 
1462  inline
1463  void complete( localIndex const iconn,
1464  StackVariables & stack ) const
1465  {
1466  GEOS_UNUSED_VAR( iconn );
1467  using namespace compositionalMultiphaseUtilities;
1468  if( stack.numConnectedElems ==1 )
1469  {
1470  // Setup Jacobian global row indicies
1471  // equations for COMPONENT + ENERGY balances
1472  globalIndex oneSidedEqnRowIndices[numEqn]{};
1473  for( integer ic = 0; ic < NC; ++ic )
1474  {
1475  oneSidedEqnRowIndices[ic] = stack.offsetUp + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1476  }
1477 
1478  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1479  globalIndex oneSidedDofColIndices_dPresCompTempUp[CP_Deriv::nDer]{};
1480  globalIndex oneSidedDofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1481  // Note localFluxJacobian cols are stored using CP_Deriv order (dP dC or dP dT dC)
1482  int ioff=0;
1483  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1484 
1485  if constexpr ( IS_THERMAL )
1486  {
1487  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1488  }
1489  for( integer jdof = 0; jdof < NC; ++jdof )
1490  {
1491  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1492  }
1493  if( m_useTotalMassEquation > 0 )
1494  {
1495  // Apply equation/variable change transformation(s)
1496  real64 work[CP_Deriv::nDer]{};
1497  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, 1, stack.localFluxJacobian_dQ, work );
1498  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, CP_Deriv::nDer, stack.localFluxJacobian, work );
1499  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, stack.localFlux );
1500  }
1501  for( integer i = 0; i < numEqn; ++i )
1502  {
1503  if( oneSidedEqnRowIndices[i] >= 0 && oneSidedEqnRowIndices[i] < m_localMatrix.numRows() )
1504  {
1505  m_localMatrix.addToRow< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1506  &oneSidedDofColIndices_dRate,
1507  stack.localFluxJacobian_dQ[i],
1508  1 );
1509  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1510  oneSidedDofColIndices_dPresCompTempUp,
1511  stack.localFluxJacobian[i],
1512  CP_Deriv::nDer );
1513  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[oneSidedEqnRowIndices[i]], stack.localFlux[i] );
1514  }
1515  }
1516  }
1517  else
1518  {
1519  // Setup Jacobian global row indicies
1520  // equations for COMPONENT + ENERGY balances
1521  globalIndex eqnRowIndices[2*numEqn]{};
1522 
1523  for( integer ic = 0; ic < NC; ++ic )
1524  {
1525  // mass balance equations for all components
1526  eqnRowIndices[TAG::NEXT *numEqn+ic] = stack.offsetNext + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1527  eqnRowIndices[TAG::CURRENT *numEqn+ic] = stack.offsetCurrent + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1528  }
1529 
1530  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1531  globalIndex dofColIndices_dPresCompUp[CP_Deriv::nDer]{};
1532  globalIndex dofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1533 
1534  int ioff=0;
1535  // Indice storage order reflects local jac col storage order CP::Deriv order P T DENS
1536  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1537 
1538  if constexpr ( IS_THERMAL )
1539  {
1540  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1541  }
1542  for( integer jdof = 0; jdof < NC; ++jdof )
1543  {
1544  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1545  }
1546 
1547 
1548  if( m_useTotalMassEquation > 0 )
1549  {
1550  // Apply equation/variable change transformation(s)
1551  real64 work[CP_Deriv::nDer]{};
1552  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, 1, 2, stack.localFluxJacobian_dQ, work );
1553  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, CP_Deriv::nDer, 2, stack.localFluxJacobian, work );
1554  shiftBlockElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, numEqn, 2, stack.localFlux );
1555  }
1556  // Note this updates diag and offdiag
1557  for( integer i = 0; i < 2*NC; ++i )
1558  {
1559  if( eqnRowIndices[i] >= 0 && eqnRowIndices[i] < m_localMatrix.numRows() )
1560  {
1561  m_localMatrix.addToRow< parallelDeviceAtomic >( eqnRowIndices[i],
1562  &dofColIndices_dRate,
1563  stack.localFluxJacobian_dQ[i],
1564  1 );
1565  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( eqnRowIndices[i],
1566  dofColIndices_dPresCompUp,
1567  stack.localFluxJacobian[i],
1568  CP_Deriv::nDer );
1569  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[eqnRowIndices[i]], stack.localFlux[i] );
1570 
1571  }
1572  }
1573  }
1574  }
1575 
1577  inline
1578  void
1579  computeExit( real64 const & dt,
1580  real64 const ( &compFlux )[NC ],
1581  StackVariables & stack,
1582  real64 ( & dCompFlux)[NC][numDof] ) const
1583  {
1584  for( integer ic = 0; ic < NC; ++ic )
1585  {
1586  stack.localFlux[ic] = -dt * compFlux[ic];
1587  // derivative with respect to rate
1588  stack.localFluxJacobian_dQ[ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1589  // derivative with respect to upstream pressure
1590  stack.localFluxJacobian[ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1591  // derivatives with respect to upstream component densities
1592  for( integer jdof = 0; jdof < NC; ++jdof )
1593  {
1594  stack.localFluxJacobian[ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1595  }
1596  if constexpr ( IS_THERMAL )
1597  {
1598  stack.localFluxJacobian[ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1599  }
1600  }
1601  }
1602 
1604  inline
1605  void
1606  compute( real64 const & dt,
1607  real64 const ( &compFlux )[NC ],
1608  StackVariables & stack,
1609  real64 ( & dCompFlux)[(NC )][numDof] ) const
1610  {
1611  // flux terms
1612  for( integer ic = 0; ic < NC; ++ic )
1613  {
1614  stack.localFlux[TAG::NEXT * NC +ic] = dt * compFlux[ic];
1615  stack.localFlux[TAG::CURRENT * NC +ic] = -dt * compFlux[ic];
1616  // derivative with respect to rate
1617  stack.localFluxJacobian_dQ[TAG::NEXT * NC+ ic][0] = dt * dCompFlux[ic][WJ_COFFSET::dQ];
1618  stack.localFluxJacobian_dQ[TAG::CURRENT * NC + +ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1619 
1620  // derivative with respect to upstream pressure
1621  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dP] = dt * dCompFlux[ic][WJ_COFFSET::dP];
1622  stack.localFluxJacobian[TAG::CURRENT * NC+ ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1623 
1624  if constexpr ( IS_THERMAL )
1625  {
1626  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dT] = dt * dCompFlux[ic][WJ_COFFSET::dT];
1627  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1628  }
1629 
1630  // derivatives with respect to upstream component densities
1631  for( integer jdof = 0; jdof < NC; ++jdof )
1632  {
1633  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dC+jdof] = dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1634  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1635  }
1636  }
1637  }
1638 
1646  template< typename FUNC = NoOpFunc >
1648  inline
1649  void computeFlux( localIndex const iwelem,
1650  StackVariables & stack,
1651  FUNC && compFluxKernelOp = NoOpFunc{} ) const
1652  {
1653 
1654  using namespace compositionalMultiphaseUtilities;
1655 
1656  // create local work arrays
1657  real64 compFracUp[NC]{};
1658  real64 compFlux[NC]{};
1659  real64 dComp[NC][NC];
1660  real64 dCompFlux[NC][numDof]{};
1661  for( integer ic = 0; ic < NC; ++ic )
1662  {
1663  for( integer jc = 0; jc < NC; ++jc )
1664  {
1665  dComp[ic][jc]=0.0;
1666  }
1667  }
1668  // Step 1) decide the upwind well element
1669 
1670  /* currentConnRate < 0 flow from iwelem to iwelemNext
1671  * currentConnRate > 0 flow from iwelemNext to iwelem
1672  * With this convention, currentConnRate < 0 at the last connection for a producer
1673  * currentConnRate > 0 at the last connection for a injector
1674  */
1675 
1676  localIndex iwelemNext = m_nextWellElemIndex[iwelem];
1677 
1678  // Current element is open and next element is closed => no dependency on next segment
1679  if( iwelemNext >= 0 )
1680  {
1681  if( m_elemStatus[iwelemNext]==WellElementSubRegion::WellElemStatus::CLOSED )
1682  {
1683  iwelemNext = -1;
1684  }
1685  }
1686 
1687  real64 const currentConnRate = m_connRate[iwelem];
1688  localIndex iwelemUp = -1;
1689  if( iwelemNext < 0 && !m_isProducer ) // exit connection, injector
1690  {
1691  // we still need to define iwelemUp for Jacobian assembly
1692  iwelemUp = iwelem;
1693 
1694  // just copy the injection stream into compFrac
1695  for( integer ic = 0; ic < NC; ++ic )
1696  {
1697  compFracUp[ic] = m_injection[ic];
1698  for( integer jc = 0; jc < NC; ++jc )
1699  {
1700  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1701  }
1702  for( integer jc = 0; jc < NC; ++jc )
1703  {
1704  dCompFlux[ic][WJ_COFFSET::dC+jc] = 0.0;
1705  }
1706  }
1707  }
1708  else
1709  {
1710  // first set iwelemUp to the upstream cell
1711  if( ( iwelemNext < 0 && m_isProducer ) // exit connection, producer
1712  || currentConnRate < 0 ) // not an exit connection, iwelem is upstream
1713  {
1714  iwelemUp = iwelem;
1715  }
1716  else // not an exit connection, iwelemNext is upstream
1717  {
1718  iwelemUp = iwelemNext;
1719  }
1720 
1721  // copy the vars of iwelemUp into compFrac
1722  for( integer ic = 0; ic < NC; ++ic )
1723  {
1724  compFracUp[ic] = m_wellElemCompFrac[iwelemUp][ic];
1725  for( integer jc = 0; jc < NC; ++jc )
1726  {
1727  dCompFlux[ic][WJ_COFFSET::dC+jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1728  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1729  }
1730  }
1731  }
1732 
1733  // Step 2) compute upstream transport coefficient
1734 
1735  for( integer ic = 0; ic < NC; ++ic )
1736  {
1737  compFlux[ic] = compFracUp[ic] * currentConnRate;
1738  dCompFlux[ic][WJ_COFFSET::dQ] = compFracUp[ic];
1739  // none of these quantities depend on pressure
1740  dCompFlux[ic][WJ_COFFSET::dP] = 0.0;
1741  if constexpr ( IS_THERMAL )
1742  {
1743  dCompFlux[ic][WJ_COFFSET::dT] = 0.0;
1744  }
1745  for( integer jc = 0; jc < NC; ++jc )
1746  {
1747  dCompFlux[ic][WJ_COFFSET::dC+jc] = dCompFlux[ic][WJ_COFFSET::dC+jc] * currentConnRate;
1748  }
1749  }
1750 
1751  stack.offsetUp = m_wellElemDofNumber[iwelemUp];
1752  stack.iwelemUp = iwelemUp;
1753  stack.offsetCurrent = m_wellElemDofNumber[iwelem];
1754  stack.iwelemCurrent= iwelem;
1755 
1756 
1757  if( iwelemNext < 0 ) // exit connection
1758  {
1759  // for this case, we only need NC mass conservation equations
1760  computeExit ( m_dt,
1761  compFlux,
1762  stack,
1763  dCompFlux );
1764 
1765  }
1766  else // not an exit connection
1767  {
1768  compute( m_dt,
1769  compFlux,
1770  stack,
1771  dCompFlux
1772  );
1773  stack.offsetNext = m_wellElemDofNumber[iwelemNext];
1774  }
1775  stack.iwelemNext = iwelemNext;
1776  compFluxKernelOp( iwelemNext, iwelemUp, currentConnRate, dComp );
1777  }
1778 
1779 
1787  template< typename POLICY, typename KERNEL_TYPE >
1788  static void
1789  launch( localIndex const numElements,
1790  KERNEL_TYPE const & kernelComponent )
1791  {
1793  forAll< POLICY >( numElements, [=] GEOS_HOST_DEVICE ( localIndex const ie )
1794  {
1795  typename KERNEL_TYPE::StackVariables stack( 1 );
1796  if( kernelComponent.skipElement( ie ))
1797  {
1798  return;
1799  }
1800  kernelComponent.setup( ie, stack );
1801  kernelComponent.computeFlux( ie, stack );
1802  kernelComponent.complete( ie, stack );
1803  } );
1804  }
1805 
1806 protected:
1808  real64 const m_dt;
1811 
1816 
1819 
1822 
1823 
1828 
1833 
1836 
1838  bool const m_isProducer;
1839 
1842 
1843 
1844 };
1845 
1850 {
1851 public:
1852 
1866  template< typename POLICY >
1867  static void
1868  createAndLaunch( integer const numComps,
1869  real64 const dt,
1870  globalIndex const rankOffset,
1871  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags,
1872  string const dofKey,
1873  WellControls const & wellControls,
1874  WellElementSubRegion const & subRegion,
1875  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1876  arrayView1d< real64 > const & localRhs )
1877  {
1878  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComps, [&]( auto NC )
1879  {
1880  integer constexpr NUM_COMP = NC();
1881 
1882  using kernelType = FaceBasedAssemblyKernel< NUM_COMP, 0 >;
1883  kernelType kernel( dt, rankOffset, dofKey, wellControls, subRegion, localMatrix, localRhs, kernelFlags );
1884  kernelType::template launch< POLICY >( subRegion.size(), kernel );
1885  } );
1886  }
1887 };
1888 } // end namespace compositionalMultiphaseWellKernels
1889 
1890 } // end namespace geos
1891 
1892 #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.
real64 getConstraintValue(real64 const &currentTime) const
Get the target bottom hole pressure value.
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.
integer getConstraintPhaseIndex() const
Const accessor for the phase constraint index.
MinimumBHPConstraint * getMinBHPConstraint()
Getters for constraints.
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, 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.