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 Total, S.A
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"
36 #include "physicsSolvers/fluidFlow/kernels/compositional/AccumulationKernel.hpp"
41 #include "physicsSolvers/fluidFlow/wells/WellControls.hpp"
43 
44 namespace geos
45 {
46 
47 
48 namespace compositionalMultiphaseWellKernels
49 {
50 
51 static constexpr real64 minDensForDivision = 1e-10;
52 
53 // tag to access well and reservoir elements in perforation rates computation
55 {
56  static constexpr integer RES = 0;
57  static constexpr integer WELL = 1;
58 };
59 
60 // tag to access the next and current well elements of a connection
61 struct ElemTag
62 {
63  static constexpr integer CURRENT = 0;
64  static constexpr integer NEXT = 1;
65 };
66 
67 // define the column offset of the derivatives
68 struct ColOffset
69 {
70  static constexpr integer DPRES = 0;
71  static constexpr integer DCOMP = 1;
72 };
73 
74 template< integer NC, integer IS_THERMAL >
76 
77 template< integer NC >
78 struct ColOffset_WellJac< NC, 0 >
79 {
80  static constexpr integer dP = 0;
81  static constexpr integer dC = 1;
82  static constexpr integer dQ = dC + NC;
83  static integer constexpr nDer = dQ + 1;
84 
85 };
86 
87 template< integer NC >
88 struct ColOffset_WellJac< NC, 1 >
89 {
90  static constexpr integer dP = 0;
91  static constexpr integer dC = 1;
92  static constexpr integer dQ = dC + NC;
93  static constexpr integer dT = dQ+1;
95  static integer constexpr nDer = dT + 1;
96 };
97 
98 // define the row offset of the residual equations
99 struct RowOffset
100 {
101  static constexpr integer CONTROL = 0;
102  static constexpr integer MASSBAL = 1;
103 };
104 
105 template< integer NC, integer IS_THERMAL >
107 
108 template< integer NC >
109 struct RowOffset_WellJac< NC, 0 >
110 {
111  static constexpr integer CONTROL = 0;
112  static constexpr integer MASSBAL = 1;
113  static constexpr integer VOLBAL = MASSBAL + NC;
114  static constexpr integer nEqn = VOLBAL+1;
115 };
116 
117 template< integer NC >
118 struct RowOffset_WellJac< NC, 1 >
119 {
120  static constexpr integer CONTROL = 0;
121  static constexpr integer MASSBAL = 1;
122  static constexpr integer VOLBAL = MASSBAL + NC;
123  static constexpr integer ENERGYBAL = VOLBAL+1;
124  static constexpr integer nEqn = ENERGYBAL+1;
125 
126 };
127 /******************************** ControlEquationHelper ********************************/
129 {
132 
134  inline
135  static
136  void
137  switchControl( bool const isProducer,
138  WellControls::Control const & 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  WellControls::Control & newControl );
148 
149  template< integer NC, integer IS_THERMAL >
151  inline
152  static void
153  compute( globalIndex const rankOffset,
154  WellControls::Control const currentControl,
155  integer const targetPhaseIndex,
156  real64 const & targetBHP,
157  real64 const & targetPhaseRate,
158  real64 const & targetTotalRate,
159  real64 const & targetMassRate,
160  real64 const & currentBHP,
161  arrayView1d< real64 const > const & dCurrentBHP,
162  arrayView1d< real64 const > const & currentPhaseVolRate,
163  arrayView2d< real64 const > const & dCurrentPhaseVolRate,
164  real64 const & currentTotalVolRate,
165  arrayView1d< real64 const > const & dCurrentTotalVolRate,
166  real64 const & massDensity,
167  globalIndex const dofNumber,
168  CRSMatrixView< real64, globalIndex const > const & localMatrix,
169  arrayView1d< real64 > const & localRhs );
170 
171 };
172 
173 /******************************** PressureRelationKernel ********************************/
174 
176 {
177  using Deriv = constitutive::multifluid::DerivativeOffset;
181 
182  template< integer NC, integer IS_THERMAL >
184  inline
185  static void
186  compute( real64 const & gravCoef,
187  real64 const & gravCoefNext,
188  real64 const & pres,
189  real64 const & presNext,
190  real64 const & totalMassDens,
191  real64 const & totalMassDensNext,
194  real64 & localPresRel,
195  real64 ( &localPresRelJacobian )[2*(NC+1+IS_THERMAL)] );
196 
197  template< integer NC, integer IS_THERMAL >
198  static void
199  launch( localIndex const size,
200  globalIndex const rankOffset,
201  bool const isLocallyOwned,
202  localIndex const iwelemControl,
203  integer const targetPhaseIndex,
204  WellControls const & wellControls,
205  real64 const & timeAtEndOfStep,
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  localIndex const numPerforations,
284  WellControls const & wellControls,
285  real64 const & currentTime,
291  arrayView1d< localIndex const > const & resElementRegion,
292  arrayView1d< localIndex const > const & resElementSubRegion,
293  arrayView1d< localIndex const > const & resElementIndex,
294  arrayView1d< real64 const > const & perfGravCoef,
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  integer const targetPhaseIndex,
324  WellControls const & wellControls,
325  real64 const & currentTime,
328  arrayView1d< real64 > const & connRate );
329 
330 };
331 
332 
333 /******************************** TotalMassDensityKernel ****************************/
334 
341 template< integer NUM_COMP, integer NUM_PHASE >
343 {
344 public:
345 
347  using Base::numComp;
348 
350  static constexpr integer numPhase = NUM_PHASE;
351 
358  constitutive::MultiFluidBase const & fluid )
359  : Base(),
360  m_phaseVolFrac( subRegion.getField< fields::well::phaseVolumeFraction >() ),
361  m_dPhaseVolFrac( subRegion.getField< fields::well::dPhaseVolumeFraction >() ),
362  m_dCompFrac_dCompDens( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
363  m_phaseMassDens( fluid.phaseMassDensity() ),
364  m_dPhaseMassDens( fluid.dPhaseMassDensity() ),
365  m_totalMassDens( subRegion.getField< fields::well::totalMassDensity >() ),
366  m_dTotalMassDens( subRegion.getField< fields::well::dTotalMassDensity >() )
367  {}
368 
375  template< typename FUNC = NoOpFunc >
377  inline
378  void compute( localIndex const ei,
379  FUNC && totalMassDensityKernelOp = NoOpFunc{} ) const
380  {
381  using Deriv = constitutive::multifluid::DerivativeOffset;
382 
383  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
384  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
385  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
386  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseMassDens = m_phaseMassDens[ei][0];
387  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseMassDens = m_dPhaseMassDens[ei][0];
388  real64 & totalMassDens = m_totalMassDens[ei];
389  arraySlice1d< real64, compflow::USD_FLUID_DC - 1 > dTotalMassDens = m_dTotalMassDens[ei];
390 
391  real64 dMassDens_dC[numComp]{};
392 
393  totalMassDens = 0.0;
394 
395  dTotalMassDens[Deriv::dP]=0.0;
396  for( integer ic = 0; ic < numComp; ++ic )
397  {
398  dTotalMassDens[Deriv::dC+ic]=0.0;
399  }
400 
401  for( integer ip = 0; ip < numPhase; ++ip )
402  {
403  totalMassDens += phaseVolFrac[ip] * phaseMassDens[ip];
404  dTotalMassDens[Deriv::dP] += dPhaseVolFrac[ip][Deriv::dP] * phaseMassDens[ip] + phaseVolFrac[ip] * dPhaseMassDens[ip][Deriv::dP];
405 
406  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseMassDens[ip], dMassDens_dC, Deriv::dC );
407  for( integer ic = 0; ic < numComp; ++ic )
408  {
409  dTotalMassDens[Deriv::dC+ic] += dPhaseVolFrac[ip][Deriv::dC+ic] * phaseMassDens[ip]
410  + phaseVolFrac[ip] * dMassDens_dC[ic];
411  }
412 
413  totalMassDensityKernelOp( ip ); //, phaseVolFrac, dTotalMassDens_dPres, dTotalMassDens_dCompDens );
414  }
415 
416  }
417 
418 protected:
419 
420  // inputs
421 
426 
430 
431  // outputs
432 
436 
437 
438 };
439 
444 {
445 public:
446 
455  template< typename POLICY >
456  static void
457  createAndLaunch( integer const numComp,
458  integer const numPhase,
459  ObjectManagerBase & subRegion,
460  constitutive::MultiFluidBase const & fluid )
461  {
462  if( numPhase == 2 )
463  {
464  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
465  {
466  integer constexpr NUM_COMP = NC();
467  TotalMassDensityKernel< NUM_COMP, 2 > kernel( subRegion, fluid );
468  TotalMassDensityKernel< NUM_COMP, 2 >::template launch< POLICY >( subRegion.size(), kernel );
469  } );
470  }
471  else if( numPhase == 3 )
472  {
473  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComp, [&] ( auto NC )
474  {
475  integer constexpr NUM_COMP = NC();
476  TotalMassDensityKernel< NUM_COMP, 3 > kernel( subRegion, fluid );
477  TotalMassDensityKernel< NUM_COMP, 3 >::template launch< POLICY >( subRegion.size(), kernel );
478  } );
479  }
480  }
481 };
482 
483 
484 /******************************** ResidualNormKernel ********************************/
485 
490 {
491 public:
492 
494  using Base::m_minNormalizer;
495  using Base::m_rankOffset;
496  using Base::m_localResidual;
497  using Base::m_dofNumber;
498 
499  ResidualNormKernel( globalIndex const rankOffset,
500  arrayView1d< real64 const > const & localResidual,
501  arrayView1d< globalIndex const > const & dofNumber,
503  integer const numComp,
504  integer const numDof,
505  integer const targetPhaseIndex,
506  WellElementSubRegion const & subRegion,
507  constitutive::MultiFluidBase const & fluid,
508  WellControls const & wellControls,
509  real64 const timeAtEndOfStep,
510  real64 const dt,
511  real64 const minNormalizer )
512  : Base( rankOffset,
513  localResidual,
514  dofNumber,
515  ghostRank,
516  minNormalizer ),
517  m_numComp( numComp ),
518  m_numDof( numDof ),
519  m_targetPhaseIndex( targetPhaseIndex ),
520  m_dt( dt ),
521  m_isLocallyOwned( subRegion.isLocallyOwned() ),
523  m_isProducer( wellControls.isProducer() ),
524  m_currentControl( wellControls.getControl() ),
525  m_targetBHP( wellControls.getTargetBHP( timeAtEndOfStep ) ),
526  m_targetTotalRate( wellControls.getTargetTotalRate( timeAtEndOfStep ) ),
527  m_targetPhaseRate( wellControls.getTargetPhaseRate( timeAtEndOfStep ) ),
528  m_targetMassRate( wellControls.getTargetMassRate( timeAtEndOfStep ) ),
529  m_volume( subRegion.getElementVolume() ),
530  m_phaseDens_n( fluid.phaseDensity_n() ),
531  m_totalDens_n( fluid.totalDensity_n() )
532  {}
533 
535  virtual void computeLinf( localIndex const iwelem,
536  LinfStackVariables & stack ) const override
537  {
539 
540  real64 normalizer = 0.0;
541  for( integer idof = 0; idof < m_numDof; ++idof )
542  {
543 
544  // Step 1: compute a normalizer for the control or pressure equation
545 
546  // for the control equation, we distinguish two cases
547  if( idof == ROFFSET::CONTROL )
548  {
549 
550  // for the top well element, normalize using the current control
551  if( m_isLocallyOwned && iwelem == m_iwelemControl )
552  {
554  {
555  // the residual entry is in pressure units
556  normalizer = m_targetBHP;
557  }
559  {
560  // the residual entry is in volume / time units
561  normalizer = LvArray::math::max( LvArray::math::abs( m_targetTotalRate ), m_minNormalizer );
562  }
564  {
565  // the residual entry is in volume / time units
566  normalizer = LvArray::math::max( LvArray::math::abs( m_targetPhaseRate ), m_minNormalizer );
567  }
569  {
570  // the residual entry is in volume / time units
571  normalizer = LvArray::math::max( LvArray::math::abs( m_targetMassRate ), m_minNormalizer );
572  }
573  }
574  // for the pressure difference equation, always normalize by the BHP
575  else
576  {
577  normalizer = m_targetBHP;
578  }
579  }
580  // Step 2: compute a normalizer for the mass balance equations
581  else if( idof >= ROFFSET::MASSBAL && idof < ROFFSET::MASSBAL + m_numComp )
582  {
583  if( m_isProducer ) // only PHASEVOLRATE is supported for now
584  {
585  // the residual is in mass units
586  normalizer = m_dt * LvArray::math::abs( m_targetPhaseRate ) * m_phaseDens_n[iwelem][0][m_targetPhaseIndex];
587  }
588  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
589  {
591  {
592  normalizer = m_dt * LvArray::math::abs( m_targetMassRate );
593  }
594  else
595  {
596  // the residual is in mass units
597  normalizer = m_dt * LvArray::math::abs( m_targetTotalRate ) * m_totalDens_n[iwelem][0];
598  }
599 
600  }
601 
602  // to make sure that everything still works well if the rate is zero, we add this check
603  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] * m_totalDens_n[iwelem][0] );
604  }
605  // Step 3: compute a normalizer for the volume balance equations
606  else if( idof == ROFFSET::MASSBAL + m_numComp )
607  {
608  if( m_isProducer ) // only PHASEVOLRATE is supported for now
609  {
610  // the residual is in volume units
611  normalizer = m_dt * LvArray::math::abs( m_targetPhaseRate );
612  }
613  else // Type::INJECTOR, only TOTALVOLRATE is supported for now
614  {
616  {
617  normalizer = m_dt * LvArray::math::abs( m_targetMassRate/ m_totalDens_n[iwelem][0] );
618  }
619  else
620  {
621  normalizer = m_dt * LvArray::math::abs( m_targetTotalRate );
622  }
623 
624  }
625 
626  }
627 
628  // to make sure that everything still works well if the rate is zero, we add this check
629  normalizer = LvArray::math::max( normalizer, m_volume[iwelem] );
630 
631  // Step 4: compute the contribution to the residual
632  real64 const val = LvArray::math::abs( m_localResidual[stack.localRow + idof] ) / normalizer;
633  if( val > stack.localValue[0] )
634  {
635  stack.localValue[0] = val;
636  }
637  }
638  }
639 
641  virtual void computeL2( localIndex const iwelem,
642  L2StackVariables & stack ) const override
643  {
644  GEOS_UNUSED_VAR( iwelem, stack );
645  GEOS_ERROR( "The L2 norm is not implemented for CompositionalMultiphaseWell" );
646  }
647 
648 
649 protected:
650 
653 
656 
659 
661  real64 const m_dt;
662 
664  bool const m_isLocallyOwned;
665 
668 
670  bool const m_isProducer;
671 
674  real64 const m_targetBHP;
675  real64 const m_targetTotalRate;
676  real64 const m_targetPhaseRate;
677  real64 const m_targetMassRate;
678 
681 
685 
686 };
687 
692 {
693 public:
694 
711  template< typename POLICY >
712  static void
713  createAndLaunch( integer const numComp,
714  integer const numDof,
715  integer const targetPhaseIndex,
716  globalIndex const rankOffset,
717  string const & dofKey,
718  arrayView1d< real64 const > const & localResidual,
719  WellElementSubRegion const & subRegion,
720  constitutive::MultiFluidBase const & fluid,
721  WellControls const & wellControls,
722  real64 const timeAtEndOfStep,
723  real64 const dt,
724  real64 const minNormalizer,
725  real64 (& residualNorm)[1] )
726  {
727  arrayView1d< globalIndex const > const dofNumber = subRegion.getReference< array1d< globalIndex > >( dofKey );
728  arrayView1d< integer const > const ghostRank = subRegion.ghostRank();
729 
730  ResidualNormKernel kernel( rankOffset, localResidual, dofNumber, ghostRank,
731  numComp, numDof, targetPhaseIndex, subRegion, fluid, wellControls, timeAtEndOfStep, dt, minNormalizer );
732  ResidualNormKernel::launchLinf< POLICY >( subRegion.size(), kernel, residualNorm );
733  }
734 
735 };
736 
737 /******************************** SolutionScalingKernel ********************************/
738 
743 {
744 public:
745 
759  template< typename POLICY >
761  createAndLaunch( real64 const maxRelativePresChange,
762  real64 const maxAbsolutePresChange,
763  real64 const maxCompFracChange,
764  real64 const maxRelativeCompDensChange,
765  globalIndex const rankOffset,
766  integer const numComp,
767  string const dofKey,
768  ElementSubRegionBase & subRegion,
769  arrayView1d< real64 const > const localSolution )
770  {
771  arrayView1d< real64 const > const pressure =
772  subRegion.getField< fields::well::pressure >();
774  subRegion.getField< fields::well::globalCompDensity >();
775  arrayView1d< real64 > pressureScalingFactor =
776  subRegion.getField< fields::well::pressureScalingFactor >();
777  arrayView1d< real64 > compDensScalingFactor =
778  subRegion.getField< fields::well::globalCompDensityScalingFactor >();
779  isothermalCompositionalMultiphaseBaseKernels::
780  SolutionScalingKernel kernel( maxRelativePresChange, maxAbsolutePresChange, maxCompFracChange, maxRelativeCompDensChange, rankOffset,
781  numComp, dofKey, subRegion, localSolution, pressure, compDens, pressureScalingFactor, compDensScalingFactor );
782  return isothermalCompositionalMultiphaseBaseKernels::
783  SolutionScalingKernel::
784  launch< POLICY >( subRegion.size(), kernel );
785  }
786 
787 };
788 
789 /******************************** SolutionCheckKernel ********************************/
790 
795 {
796 public:
797 
809  template< typename POLICY >
811  createAndLaunch( integer const allowCompDensChopping,
812  CompositionalMultiphaseFVM::ScalingType const scalingType,
813  real64 const scalingFactor,
814  arrayView1d< real64 const > const pressure,
816  arrayView1d< real64 > pressureScalingFactor,
817  arrayView1d< real64 > compDensScalingFactor,
818  globalIndex const rankOffset,
819  integer const numComp,
820  string const dofKey,
821  ElementSubRegionBase & subRegion,
822  arrayView1d< real64 const > const localSolution )
823  {
824 
825  isothermalCompositionalMultiphaseBaseKernels::
826  SolutionCheckKernel kernel( allowCompDensChopping, 0, scalingType, scalingFactor,
827  pressure, compDens, pressureScalingFactor, compDensScalingFactor, rankOffset,
828  numComp, dofKey, subRegion, localSolution );
829  return isothermalCompositionalMultiphaseBaseKernels::
830  SolutionCheckKernel::
831  launch< POLICY >( subRegion.size(), kernel );
832  }
833 
834 };
835 
836 /******************************** ElementBasedAssemblyKernel ********************************/
837 
844 template< integer NUM_COMP, integer IS_THERMAL >
846 {
847 public:
850 
851  // Well jacobian column and row indicies
852  // tjb - change NUM_DOF to IS_THERMAL
853  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NUM_COMP, IS_THERMAL >;
857  static constexpr integer numComp = NUM_COMP;
858 
860  static constexpr integer numDof = NUM_COMP + 1 + IS_THERMAL;
861 
863  static constexpr integer numEqn = NUM_COMP + 1 + IS_THERMAL;
864 
865 
878  integer const isProducer,
879  globalIndex const rankOffset,
880  string const dofKey,
881  WellElementSubRegion const & subRegion,
882  constitutive::MultiFluidBase const & fluid,
883  CRSMatrixView< real64, globalIndex const > const & localMatrix,
884  arrayView1d< real64 > const & localRhs,
885  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const kernelFlags )
886  : m_numPhases( numPhases ),
887  m_isProducer( isProducer ),
888  m_rankOffset( rankOffset ),
889  m_iwelemControl( subRegion.getTopWellElementIndex() ),
890  m_dofNumber( subRegion.getReference< array1d< globalIndex > >( dofKey ) ),
891  m_elemGhostRank( subRegion.ghostRank() ),
892  m_volume( subRegion.getElementVolume() ),
893  m_dCompFrac_dCompDens( subRegion.getField< fields::flow::dGlobalCompFraction_dGlobalCompDensity >() ),
894  m_phaseVolFrac_n( subRegion.getField< fields::flow::phaseVolumeFraction_n >() ),
895  m_phaseVolFrac( subRegion.getField< fields::flow::phaseVolumeFraction >() ),
896  m_dPhaseVolFrac( subRegion.getField< fields::flow::dPhaseVolumeFraction >() ),
897  m_phaseDens_n( fluid.phaseDensity_n() ),
898  m_phaseDens( fluid.phaseDensity() ),
899  m_dPhaseDens( fluid.dPhaseDensity() ),
900  m_phaseCompFrac_n( fluid.phaseCompFraction_n() ),
901  m_phaseCompFrac( fluid.phaseCompFraction() ),
902  m_dPhaseCompFrac( fluid.dPhaseCompFraction() ),
903  m_compDens( subRegion.getField< fields::flow::globalCompDensity >() ),
904  m_compDens_n( subRegion.getField< fields::flow::globalCompDensity_n >() ),
905  m_localMatrix( localMatrix ),
906  m_localRhs( localRhs ),
907  m_kernelFlags( kernelFlags )
908  {}
909 
915  {
916 public:
917 
918  // volume information (used by both accumulation and volume balance)
919  real64 volume = 0.0;
920 
921  // Residual information
922 
925 
927  globalIndex dofIndices[numDof]{}; // NC compdens + P + thermal
928  globalIndex eqnRowIndices[numDof]{};
929  globalIndex dofColIndices[numDof]{};
930 
933 
936 
937  };
938 
945  integer elemGhostRank( localIndex const ei ) const
946  { return m_elemGhostRank( ei ); }
947 
948 
955  void setup( localIndex const ei,
956  StackVariables & stack ) const
957  {
958  // initialize the volume
959  stack.volume = m_volume[ei];
960 
961  // Note row/col indices needed to be consistent with layout of stack.localJacobian
962  // Setup row equation indices for this element ( mass + vol + thermal if valid)
963 
964  // 1) Mass Balance
965  for( integer ic = 0; ic < numComp; ++ic )
966  {
967  stack.eqnRowIndices[ic] = m_dofNumber[ei] + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
968  }
969 // 2) Volume Balance
970  stack.eqnRowIndices[numComp] = m_dofNumber[ei] + WJ_ROFFSET::VOLBAL - m_rankOffset;
971  // 3) Energy Balance
972  if constexpr ( IS_THERMAL )
973  {
974  stack.eqnRowIndices[numComp+1] = m_dofNumber[ei] + WJ_ROFFSET::ENERGYBAL - m_rankOffset;
975  }
976  // Setup equation column indices for this element ( P + COMPDENS + THERMAL if valid)
977  stack.dofColIndices[0] = m_dofNumber[ei] + WJ_COFFSET::dP;
978  for( integer ic = 0; ic < numComp; ++ic )
979  {
980  stack.dofColIndices[ic+1] = m_dofNumber[ei] + WJ_COFFSET::dC+ic;
981  }
982  if constexpr ( IS_THERMAL )
983  {
984  stack.dofColIndices[numComp+1] = m_dofNumber[ei] + WJ_COFFSET::dT;
985  }
986  if( 1 )
987  for( integer jc = 0; jc < numEqn; ++jc )
988  {
989  stack.localResidual[jc] = 0.0;
990  for( integer ic = 0; ic < numDof; ++ic )
991  {
992  stack.localJacobian[jc][ic] = 0.0;
993  }
994 
995  }
996 
997  }
998 
1006  template< typename FUNC = NoOpFunc >
1009  StackVariables & stack,
1010  FUNC && phaseAmountKernelOp = NoOpFunc{} ) const
1011  {
1012 
1013  using Deriv = constitutive::multifluid::DerivativeOffset;
1014 
1015  // construct the slices for variables accessed multiple times
1016  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 > dCompFrac_dCompDens = m_dCompFrac_dCompDens[ei];
1017 
1018  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac_n = m_phaseVolFrac_n[ei];
1019  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
1020  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
1021 
1022  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens_n = m_phaseDens_n[ei][0];
1023  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 > phaseDens = m_phaseDens[ei][0];
1024  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 > dPhaseDens = m_dPhaseDens[ei][0];
1025 
1026  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac_n = m_phaseCompFrac_n[ei][0];
1027  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 > phaseCompFrac = m_phaseCompFrac[ei][0];
1028  arraySlice3d< real64 const, constitutive::multifluid::USD_PHASE_COMP_DC - 2 > dPhaseCompFrac = m_dPhaseCompFrac[ei][0];
1029 
1030  // temporary work arrays
1031  real64 dPhaseAmount[FLUID_PROP_COFFSET::nDer]{};
1032  real64 dPhaseAmount_dC[numComp]{};
1033  real64 dPhaseCompFrac_dC[numComp]{};
1034 
1035  // sum contributions to component accumulation from each phase
1036  for( integer ip = 0; ip < m_numPhases; ++ip )
1037  {
1038  real64 const phaseAmount = stack.volume * phaseVolFrac[ip] * phaseDens[ip];
1039  real64 const phaseAmount_n = stack.volume * phaseVolFrac_n[ip] * phaseDens_n[ip];
1040  //remove tjb
1041  real64 const dPhaseAmount_dP = stack.volume * ( dPhaseVolFrac[ip][Deriv::dP] * phaseDens[ip]
1042  + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dP] );
1043  dPhaseAmount[FLUID_PROP_COFFSET::dP]=stack.volume * ( dPhaseVolFrac[ip][Deriv::dP] * phaseDens[ip]
1044  + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dP] );
1045 
1046  // assemble density dependence
1047  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], dPhaseAmount_dC, Deriv::dC );
1048  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseDens[ip], &dPhaseAmount[FLUID_PROP_COFFSET::dC], Deriv::dC );
1049  for( integer jc = 0; jc < numComp; ++jc )
1050  {
1051  dPhaseAmount_dC[jc] = dPhaseAmount_dC[jc] * phaseVolFrac[ip]
1052  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1053  dPhaseAmount_dC[jc] *= stack.volume;
1054  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] = dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] * phaseVolFrac[ip]
1055  + phaseDens[ip] * dPhaseVolFrac[ip][Deriv::dC+jc];
1056  dPhaseAmount[FLUID_PROP_COFFSET::dC+jc] *= stack.volume;
1057  }
1058 // tjb- remove when safe
1059  for( integer ic = 0; ic < numComp; ic++ )
1060  {
1061  assert( fabs( dPhaseAmount[FLUID_PROP_COFFSET::dC+ic] -dPhaseAmount_dC[ic] ) < FLT_EPSILON );
1062 
1063  }
1064  // ic - index of component whose conservation equation is assembled
1065  // (i.e. row number in local matrix)
1066  for( integer ic = 0; ic < numComp; ++ic )
1067  {
1068  real64 const phaseCompAmount = phaseAmount * phaseCompFrac[ip][ic];
1069  real64 const phaseCompAmount_n = phaseAmount_n * phaseCompFrac_n[ip][ic];
1070 
1071  real64 const dPhaseCompAmount_dP = dPhaseAmount_dP * phaseCompFrac[ip][ic]
1072  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dP];
1073 
1074  stack.localResidual[ic] += phaseCompAmount - phaseCompAmount_n;
1075  stack.localJacobian[ic][0] += dPhaseCompAmount_dP;
1076 
1077  // jc - index of component w.r.t. whose compositional var the derivative is being taken
1078  // (i.e. col number in local matrix)
1079 
1080  // assemble phase composition dependence
1081  applyChainRule( numComp, dCompFrac_dCompDens, dPhaseCompFrac[ip][ic], dPhaseCompFrac_dC, Deriv::dC );
1082  for( integer jc = 0; jc < numComp; ++jc )
1083  {
1084  real64 const dPhaseCompAmount_dC = dPhaseCompFrac_dC[jc] * phaseAmount
1085  + phaseCompFrac[ip][ic] * dPhaseAmount[FLUID_PROP_COFFSET::dC+jc];
1086 
1087  stack.localJacobian[ic][jc + 1] += dPhaseCompAmount_dC;
1088  }
1089  }
1090  if constexpr ( IS_THERMAL )
1091  {
1092  dPhaseAmount[FLUID_PROP_COFFSET::dT] = stack.volume * (dPhaseVolFrac[ip][Deriv::dT] * phaseDens[ip] + phaseVolFrac[ip] * dPhaseDens[ip][Deriv::dT] );
1093  for( integer ic = 0; ic < numComp; ++ic )
1094  {
1095  // assemble the derivatives of the component mass balance equations with respect to temperature
1096  stack.localJacobian[ic][numComp+1] += dPhaseAmount[FLUID_PROP_COFFSET::dT] * phaseCompFrac[ip][ic]
1097  + phaseAmount * dPhaseCompFrac[ip][ic][Deriv::dT];
1098  }
1099  }
1100  // call the lambda in the phase loop to allow the reuse of the phase amounts and their derivatives
1101  // possible use: assemble accumulation term of the energy equation for this phase
1102  phaseAmountKernelOp( ip, phaseAmount, phaseAmount_n, dPhaseAmount );
1103 
1104  }
1105 
1106  // check zero diagonal (works only in debug)
1107  /*
1108  for( integer ic = 0; ic < numComp; ++ic )
1109  {
1110  GEOS_ASSERT_MSG ( LvArray::math::abs( stack.localJacobian[ic][ic] ) > minDensForDivision,
1111  GEOS_FMT( "Zero diagonal in Jacobian: equation {}, value = {}", ic, stack.localJacobian[ic][ic] ) );
1112  }
1113  */
1114  }
1115 
1116 
1127  StackVariables & stack ) const
1128  {
1129  using Deriv = constitutive::multifluid::DerivativeOffset;
1130 
1131  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > phaseVolFrac = m_phaseVolFrac[ei];
1132  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > dPhaseVolFrac = m_dPhaseVolFrac[ei];
1133 
1134  real64 oneMinusPhaseVolFracSum = 1.0;
1135 
1136  // sum contributions to component accumulation from each phase
1137 // Note localJacobian stores equation balances in order of component/vol/enerqy
1138  // These are mapped to solver orderings with indicies setup in stack variables
1139  for( integer ip = 0; ip < m_numPhases; ++ip )
1140  {
1141  oneMinusPhaseVolFracSum -= phaseVolFrac[ip];
1142  stack.localJacobian[numComp][0] -= dPhaseVolFrac[ip][Deriv::dP];
1143 
1144  for( integer jc = 0; jc < numComp; ++jc )
1145  {
1146  stack.localJacobian[numComp][jc+1] -= dPhaseVolFrac[ip][Deriv::dC+jc];
1147  }
1148 
1149  if constexpr ( IS_THERMAL)
1150  {
1151  stack.localJacobian[numComp][numComp+1] -= dPhaseVolFrac[ip][Deriv::dT];
1152  }
1153 
1154  }
1155  // scale saturation-based volume balance by pore volume (for better scaling w.r.t. other equations)
1156  stack.localResidual[numComp] = stack.volume * oneMinusPhaseVolFracSum;
1157  for( integer idof = 0; idof < numComp+1+IS_THERMAL; ++idof )
1158  {
1159  stack.localJacobian[numComp][idof] *= stack.volume;
1160  }
1161 
1162  }
1163 
1170  void complete( localIndex const ei, //GEOS_UNUSED_PARAM( ei ),
1171  StackVariables & stack ) const
1172  {
1173  using namespace compositionalMultiphaseUtilities;
1174 
1175  integer const numRows = numComp+1+ IS_THERMAL;
1176 
1177  if constexpr ( IS_THERMAL)
1178  {
1179  if( ei == m_iwelemControl && !m_isProducer )
1180  {
1181  // For top segment energy balance eqn replaced with T(n+1) - T = 0
1182  // No other energy balance derivatives
1183  // Assumption is global index == 0 is top segment with fixed temp BC
1184 
1185  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1186  {
1187  stack.localJacobian[numRows-1][i] = 0.0;
1188  }
1189  // constant Temperature
1190  for( integer i=0; i < numComp+1+IS_THERMAL; i++ )
1191  stack.localJacobian[i][numRows-1] = 0.0;
1192  stack.localJacobian[numRows-1][numRows-1] = 1.0;
1193 
1194  stack.localResidual[numRows-1]=0.0;
1195  }
1196  }
1197 
1198  if( m_kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) )
1199  {
1200  // apply equation/variable change transformation to the component mass balance equations
1201  real64 work[numComp + 1 + IS_THERMAL]{};
1202  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numComp, numComp+1+ IS_THERMAL, stack.localJacobian, work );
1203  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numComp, stack.localResidual );
1204  }
1205 
1206  // add contribution to residual and jacobian into:
1207  // - the component mass balance equations (i = 0 to i = numComp-1)
1208  // - the volume balance equations (i = numComp)
1209  // note that numDof includes derivatives wrt temperature if this class is derived in ThermalKernels
1210 
1211  for( integer i = 0; i < numRows; ++i )
1212  {
1213  m_localRhs[stack.eqnRowIndices[i]] += stack.localResidual[i];
1214  m_localMatrix.addToRow< serialAtomic >( stack.eqnRowIndices[i],
1215  stack.dofColIndices,
1216  stack.localJacobian[i],
1217  numComp+1+ IS_THERMAL );
1218  }
1219 
1220  }
1221 
1229  template< typename POLICY, typename KERNEL_TYPE >
1230  static void
1231  launch( localIndex const numElems,
1232  KERNEL_TYPE const & kernelComponent )
1233  {
1235 
1236  forAll< POLICY >( numElems, [=] GEOS_HOST_DEVICE ( localIndex const ei )
1237  {
1238  if( kernelComponent.elemGhostRank( ei ) >= 0 )
1239  {
1240  return;
1241  }
1242 
1243  typename KERNEL_TYPE::StackVariables stack;
1244 
1245  kernelComponent.setup( ei, stack );
1246  kernelComponent.computeAccumulation( ei, stack );
1247  kernelComponent.computeVolumeBalance( ei, stack );
1248  kernelComponent.complete( ei, stack );
1249  } );
1250  }
1251 
1252 protected:
1253 
1256 
1259 
1262 
1265 
1268 
1271 
1274 
1277  arrayView2d< real64 const > const m_porosity;
1278  arrayView2d< real64 const > const m_dPoro_dPres;
1279 
1282 
1287 
1292 
1297 
1298  // Views on component densities
1301 
1306 
1307  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > const m_kernelFlags;
1308 };
1309 
1310 
1315 {
1316 public:
1329  template< typename POLICY >
1330  static void
1331  createAndLaunch( localIndex const numComps,
1332  localIndex const numPhases,
1333  integer const isProducer,
1334  globalIndex const rankOffset,
1335  integer const useTotalMassEquation,
1336  string const dofKey,
1337  WellElementSubRegion const & subRegion,
1338  constitutive::MultiFluidBase const & fluid,
1339  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1340  arrayView1d< real64 > const & localRhs )
1341  {
1342  geos::internal::kernelLaunchSelectorCompThermSwitch( numComps, 0, [&]( auto NC, auto IS_THERMAL )
1343  {
1344  localIndex constexpr NUM_COMP = NC();
1345 
1346  integer constexpr istherm = IS_THERMAL();
1347 
1348  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags;
1349  if( useTotalMassEquation )
1350  kernelFlags.set( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation );
1351 
1353  kernel( numPhases, isProducer, rankOffset, dofKey, subRegion, fluid, localMatrix, localRhs, kernelFlags );
1355  launch< POLICY, ElementBasedAssemblyKernel< NUM_COMP, istherm > >( subRegion.size(), kernel );
1356  } );
1357  }
1358 };
1367 template< integer NC, integer IS_THERMAL >
1369 {
1370 public:
1371 
1375 
1376  using FLUID_PROP_COFFSET = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1379 
1380  using CP_Deriv = constitutive::multifluid::DerivativeOffsetC< NC, IS_THERMAL >;
1382  static constexpr integer numComp = NC;
1383 
1385  static constexpr integer numDof = WJ_COFFSET::nDer;
1386 
1388  static constexpr integer numEqn = NC;
1389 
1390  static constexpr integer maxNumElems = 2;
1391  static constexpr integer maxStencilSize = 2;
1407  globalIndex const rankOffset,
1408  string const wellDofKey,
1409  WellControls const & wellControls,
1410  WellElementSubRegion const & subRegion,
1411  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1412  arrayView1d< real64 > const & localRhs,
1413  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags )
1414  :
1415  m_dt( dt ),
1416  m_rankOffset( rankOffset ),
1417  m_wellElemDofNumber ( subRegion.getReference< array1d< globalIndex > >( wellDofKey ) ),
1418  m_nextWellElemIndex ( subRegion.getReference< array1d< localIndex > >( WellElementSubRegion::viewKeyStruct::nextWellElementIndexString()) ),
1419  m_connRate ( subRegion.getField< fields::well::mixtureConnectionRate >() ),
1420  m_wellElemCompFrac ( subRegion.getField< fields::well::globalCompFraction >() ),
1421  m_dWellElemCompFrac_dCompDens ( subRegion.getField< fields::well::dGlobalCompFraction_dGlobalCompDensity >() ),
1422  m_localMatrix( localMatrix ),
1423  m_localRhs ( localRhs ),
1424  m_useTotalMassEquation ( kernelFlags.isSet( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation ) ),
1425  m_isProducer ( wellControls.isProducer() ),
1426  m_injection ( wellControls.getInjectionStream() )
1427  {}
1428 
1430  {
1431 public:
1432 
1440  : stencilSize( size ),
1441  numConnectedElems( 2 ),
1442  dofColIndices( size * numDof )
1443  {}
1444 
1445  // Stencil information
1446  localIndex const stencilSize;
1449 
1450 
1451  // edge indexes
1452  localIndex iwelemUp;
1453  localIndex iwelemNext;
1454  localIndex iwelemCurrent;
1455  globalIndex offsetUp;
1456  globalIndex offsetCurrent;
1457  globalIndex offsetNext;
1458  // Local degrees of freedom and local residual/jacobian
1459 
1462 
1469  };
1470 
1477  inline
1478  void setup( localIndex const iconn,
1479  StackVariables & stack ) const
1480  {
1481  stack.numConnectedElems=2;
1482  if( m_nextWellElemIndex[iconn] <0 )
1483  {
1484  stack.numConnectedElems = 1;
1485  }
1486  stack.localFlux.resize( stack.numConnectedElems*numEqn );
1487  stack.localFluxJacobian.resize( stack.numConnectedElems * numEqn, stack.stencilSize * numDof );
1488  stack.localFluxJacobian_dQ.resize( stack.numConnectedElems * numEqn, 1 );
1489 
1490  }
1491 
1498  inline
1499  void complete( localIndex const iconn,
1500  StackVariables & stack ) const
1501  {
1502  GEOS_UNUSED_VAR( iconn );
1503  using namespace compositionalMultiphaseUtilities;
1504  if( stack.numConnectedElems ==1 )
1505  {
1506  // Setup Jacobian global row indicies
1507  // equations for COMPONENT + ENERGY balances
1508  globalIndex oneSidedEqnRowIndices[numEqn]{};
1509  for( integer ic = 0; ic < NC; ++ic )
1510  {
1511  oneSidedEqnRowIndices[ic] = stack.offsetUp + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1512  }
1513 
1514  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1515  globalIndex oneSidedDofColIndices_dPresCompTempUp[CP_Deriv::nDer]{};
1516  globalIndex oneSidedDofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1517  // Note localFluxJacobian cols are stored using CP_Deriv order (dP dC or dP dT dC)
1518  int ioff=0;
1519  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1520 
1521  if constexpr ( IS_THERMAL )
1522  {
1523  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1524  }
1525  for( integer jdof = 0; jdof < NC; ++jdof )
1526  {
1527  oneSidedDofColIndices_dPresCompTempUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1528  }
1529  if( m_useTotalMassEquation > 0 )
1530  {
1531  // Apply equation/variable change transformation(s)
1532  real64 work[CP_Deriv::nDer]{};
1533  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, 1, stack.localFluxJacobian_dQ, work );
1534  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, CP_Deriv::nDer, stack.localFluxJacobian, work );
1535  shiftElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, stack.localFlux );
1536  }
1537  for( integer i = 0; i < numEqn; ++i )
1538  {
1539  if( oneSidedEqnRowIndices[i] >= 0 && oneSidedEqnRowIndices[i] < m_localMatrix.numRows() )
1540  {
1541  m_localMatrix.addToRow< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1542  &oneSidedDofColIndices_dRate,
1543  stack.localFluxJacobian_dQ[i],
1544  1 );
1545  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( oneSidedEqnRowIndices[i],
1546  oneSidedDofColIndices_dPresCompTempUp,
1547  stack.localFluxJacobian[i],
1548  CP_Deriv::nDer );
1549  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[oneSidedEqnRowIndices[i]], stack.localFlux[i] );
1550  }
1551  }
1552  }
1553  else
1554  {
1555  // Setup Jacobian global row indicies
1556  // equations for COMPONENT + ENERGY balances
1557  globalIndex eqnRowIndices[2*numEqn]{};
1558 
1559  for( integer ic = 0; ic < NC; ++ic )
1560  {
1561  // mass balance equations for all components
1562  eqnRowIndices[TAG::NEXT *numEqn+ic] = stack.offsetNext + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1563  eqnRowIndices[TAG::CURRENT *numEqn+ic] = stack.offsetCurrent + WJ_ROFFSET::MASSBAL + ic - m_rankOffset;
1564  }
1565 
1566  // Setup Jacobian global col indicies ( Mapping from local jac order to well jac order)
1567  globalIndex dofColIndices_dPresCompUp[CP_Deriv::nDer]{};
1568  globalIndex dofColIndices_dRate = stack.offsetCurrent + WJ_COFFSET::dQ;
1569 
1570  int ioff=0;
1571  // Indice storage order reflects local jac col storage order CP::Deriv order P T DENS
1572  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dP;
1573 
1574  if constexpr ( IS_THERMAL )
1575  {
1576  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dT;
1577  }
1578  for( integer jdof = 0; jdof < NC; ++jdof )
1579  {
1580  dofColIndices_dPresCompUp[ioff++] = stack.offsetUp + WJ_COFFSET::dC+ jdof;
1581  }
1582 
1583 
1584  if( m_useTotalMassEquation > 0 )
1585  {
1586  // Apply equation/variable change transformation(s)
1587  real64 work[CP_Deriv::nDer]{};
1588  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, 1, 2, stack.localFluxJacobian_dQ, work );
1589  shiftBlockRowsAheadByOneAndReplaceFirstRowWithColumnSum( numEqn, numEqn, CP_Deriv::nDer, 2, stack.localFluxJacobian, work );
1590  shiftBlockElementsAheadByOneAndReplaceFirstElementWithSum( numEqn, numEqn, 2, stack.localFlux );
1591  }
1592  // Note this updates diag and offdiag
1593  for( integer i = 0; i < 2*NC; ++i )
1594  {
1595  if( eqnRowIndices[i] >= 0 && eqnRowIndices[i] < m_localMatrix.numRows() )
1596  {
1597  m_localMatrix.addToRow< parallelDeviceAtomic >( eqnRowIndices[i],
1598  &dofColIndices_dRate,
1599  stack.localFluxJacobian_dQ[i],
1600  1 );
1601  m_localMatrix.addToRowBinarySearchUnsorted< parallelDeviceAtomic >( eqnRowIndices[i],
1602  dofColIndices_dPresCompUp,
1603  stack.localFluxJacobian[i],
1604  CP_Deriv::nDer );
1605  RAJA::atomicAdd( parallelDeviceAtomic{}, &m_localRhs[eqnRowIndices[i]], stack.localFlux[i] );
1606 
1607  }
1608  }
1609  }
1610  }
1611 
1613  inline
1614  void
1615  computeExit( real64 const & dt,
1616  real64 const ( &compFlux )[NC ],
1617  StackVariables & stack,
1618  real64 ( & dCompFlux)[NC][numDof] ) const
1619  {
1620  for( integer ic = 0; ic < NC; ++ic )
1621  {
1622  stack.localFlux[ic] = -dt * compFlux[ic];
1623  // derivative with respect to rate
1624  stack.localFluxJacobian_dQ[ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1625  // derivative with respect to upstream pressure
1626  stack.localFluxJacobian[ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1627  // derivatives with respect to upstream component densities
1628  for( integer jdof = 0; jdof < NC; ++jdof )
1629  {
1630  stack.localFluxJacobian[ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1631  }
1632  if constexpr ( IS_THERMAL )
1633  {
1634  stack.localFluxJacobian[ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1635  }
1636  }
1637  }
1638 
1640  inline
1641  void
1642  compute( real64 const & dt,
1643  real64 const ( &compFlux )[NC ],
1644  StackVariables & stack,
1645  real64 ( & dCompFlux)[(NC )][numDof] ) const
1646  {
1647  // flux terms
1648  for( integer ic = 0; ic < NC; ++ic )
1649  {
1650  stack.localFlux[TAG::NEXT * NC +ic] = dt * compFlux[ic];
1651  stack.localFlux[TAG::CURRENT * NC +ic] = -dt * compFlux[ic];
1652  // derivative with respect to rate
1653  stack.localFluxJacobian_dQ[TAG::NEXT * NC+ ic][0] = dt * dCompFlux[ic][WJ_COFFSET::dQ];
1654  stack.localFluxJacobian_dQ[TAG::CURRENT * NC + +ic][0] = -dt * dCompFlux[ic][WJ_COFFSET::dQ];
1655 
1656  // derivative with respect to upstream pressure
1657  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dP] = dt * dCompFlux[ic][WJ_COFFSET::dP];
1658  stack.localFluxJacobian[TAG::CURRENT * NC+ ic][CP_Deriv::dP] = -dt * dCompFlux[ic][WJ_COFFSET::dP];
1659 
1660  if constexpr ( IS_THERMAL )
1661  {
1662  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dT] = dt * dCompFlux[ic][WJ_COFFSET::dT];
1663  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dT] = -dt * dCompFlux[ic][WJ_COFFSET::dT];
1664  }
1665 
1666  // derivatives with respect to upstream component densities
1667  for( integer jdof = 0; jdof < NC; ++jdof )
1668  {
1669  stack.localFluxJacobian[TAG::NEXT * NC +ic][CP_Deriv::dC+jdof] = dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1670  stack.localFluxJacobian[TAG::CURRENT * NC +ic][CP_Deriv::dC+jdof] = -dt * dCompFlux[ic][WJ_COFFSET::dC+jdof];
1671  }
1672  }
1673  }
1674 
1682  template< typename FUNC = NoOpFunc >
1684  inline
1685  void computeFlux( localIndex const iwelem,
1686  StackVariables & stack,
1687  FUNC && compFluxKernelOp = NoOpFunc{} ) const
1688  {
1689 
1690  using namespace compositionalMultiphaseUtilities;
1691 
1692  // create local work arrays
1693  real64 compFracUp[NC]{};
1694  real64 compFlux[NC]{};
1695  real64 dComp[NC][NC];
1696  real64 dCompFlux[NC][numDof]{};
1697  for( integer ic = 0; ic < NC; ++ic )
1698  {
1699  for( integer jc = 0; jc < NC; ++jc )
1700  {
1701  dComp[ic][jc]=0.0;
1702  }
1703  }
1704  // Step 1) decide the upwind well element
1705 
1706  /* currentConnRate < 0 flow from iwelem to iwelemNext
1707  * currentConnRate > 0 flow from iwelemNext to iwelem
1708  * With this convention, currentConnRate < 0 at the last connection for a producer
1709  * currentConnRate > 0 at the last connection for a injector
1710  */
1711 
1712  localIndex const iwelemNext = m_nextWellElemIndex[iwelem];
1713  real64 const currentConnRate = m_connRate[iwelem];
1714  localIndex iwelemUp = -1;
1715 
1716  if( iwelemNext < 0 && !m_isProducer ) // exit connection, injector
1717  {
1718  // we still need to define iwelemUp for Jacobian assembly
1719  iwelemUp = iwelem;
1720 
1721  // just copy the injection stream into compFrac
1722  for( integer ic = 0; ic < NC; ++ic )
1723  {
1724  compFracUp[ic] = m_injection[ic];
1725  for( integer jc = 0; jc < NC; ++jc )
1726  {
1727  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1728  }
1729  for( integer jc = 0; jc < NC; ++jc )
1730  {
1731  dCompFlux[ic][WJ_COFFSET::dC+jc] = 0.0;
1732  }
1733  }
1734  }
1735  else
1736  {
1737  // first set iwelemUp to the upstream cell
1738  if( ( iwelemNext < 0 && m_isProducer ) // exit connection, producer
1739  || currentConnRate < 0 ) // not an exit connection, iwelem is upstream
1740  {
1741  iwelemUp = iwelem;
1742  }
1743  else // not an exit connection, iwelemNext is upstream
1744  {
1745  iwelemUp = iwelemNext;
1746  }
1747 
1748  // copy the vars of iwelemUp into compFrac
1749  for( integer ic = 0; ic < NC; ++ic )
1750  {
1751  compFracUp[ic] = m_wellElemCompFrac[iwelemUp][ic];
1752  for( integer jc = 0; jc < NC; ++jc )
1753  {
1754  dCompFlux[ic][WJ_COFFSET::dC+jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1755  dComp[ic][jc] = m_dWellElemCompFrac_dCompDens[iwelemUp][ic][jc];
1756  }
1757  }
1758  }
1759 
1760  // Step 2) compute upstream transport coefficient
1761 
1762  for( integer ic = 0; ic < NC; ++ic )
1763  {
1764  compFlux[ic] = compFracUp[ic] * currentConnRate;
1765  dCompFlux[ic][WJ_COFFSET::dQ] = compFracUp[ic];
1766  // none of these quantities depend on pressure
1767  dCompFlux[ic][WJ_COFFSET::dP] = 0.0;
1768  if constexpr ( IS_THERMAL )
1769  {
1770  dCompFlux[ic][WJ_COFFSET::dT] = 0.0;
1771  }
1772  for( integer jc = 0; jc < NC; ++jc )
1773  {
1774  dCompFlux[ic][WJ_COFFSET::dC+jc] = dCompFlux[ic][WJ_COFFSET::dC+jc] * currentConnRate;
1775  }
1776  }
1777 
1778  stack.offsetUp = m_wellElemDofNumber[iwelemUp];
1779  stack.iwelemUp = iwelemUp;
1780  stack.offsetCurrent = m_wellElemDofNumber[iwelem];
1781  stack.iwelemCurrent= iwelem;
1782 
1783 
1784  if( iwelemNext < 0 ) // exit connection
1785  {
1786  // for this case, we only need NC mass conservation equations
1787  computeExit ( m_dt,
1788  compFlux,
1789  stack,
1790  dCompFlux );
1791 
1792  }
1793  else // not an exit connection
1794  {
1795  compute( m_dt,
1796  compFlux,
1797  stack,
1798  dCompFlux
1799  );
1800  stack.offsetNext = m_wellElemDofNumber[iwelemNext];
1801  }
1802  stack.iwelemNext = iwelemNext;
1803  compFluxKernelOp( iwelemNext, iwelemUp, currentConnRate, dComp );
1804  }
1805 
1806 
1814  template< typename POLICY, typename KERNEL_TYPE >
1815  static void
1816  launch( localIndex const numElements,
1817  KERNEL_TYPE const & kernelComponent )
1818  {
1820  forAll< POLICY >( numElements, [=] GEOS_HOST_DEVICE ( localIndex const ie )
1821  {
1822  typename KERNEL_TYPE::StackVariables stack( 1 );
1823 
1824  kernelComponent.setup( ie, stack );
1825  kernelComponent.computeFlux( ie, stack );
1826  kernelComponent.complete( ie, stack );
1827  } );
1828  }
1829 
1830 protected:
1832  real64 const m_dt;
1835 
1840 
1843 
1844 
1849 
1854 
1857 
1859  bool const m_isProducer;
1860 
1863 
1864 
1865 };
1866 
1871 {
1872 public:
1873 
1887  template< typename POLICY >
1888  static void
1889  createAndLaunch( integer const numComps,
1890  real64 const dt,
1891  globalIndex const rankOffset,
1892  integer const useTotalMassEquation,
1893  string const dofKey,
1894  WellControls const & wellControls,
1895  WellElementSubRegion const & subRegion,
1896  CRSMatrixView< real64, globalIndex const > const & localMatrix,
1897  arrayView1d< real64 > const & localRhs )
1898  {
1899  isothermalCompositionalMultiphaseBaseKernels::internal::kernelLaunchSelectorCompSwitch( numComps, [&]( auto NC )
1900  {
1901  integer constexpr NUM_COMP = NC();
1902 
1903 
1904  BitFlags< isothermalCompositionalMultiphaseBaseKernels::KernelFlags > kernelFlags;
1905  if( useTotalMassEquation )
1906  kernelFlags.set( isothermalCompositionalMultiphaseBaseKernels::KernelFlags::TotalMassEquation );
1907 
1908  using kernelType = FaceBasedAssemblyKernel< NUM_COMP, 0 >;
1909 
1910 
1911  kernelType kernel( dt, rankOffset, dofKey, wellControls, subRegion, localMatrix, localRhs, kernelFlags );
1912  kernelType::template launch< POLICY >( subRegion.size(), kernel );
1913  } );
1914  }
1915 };
1916 } // end namespace compositionalMultiphaseWellKernels
1917 
1918 } // end namespace geos
1919 
1920 #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, integer const useTotalMassEquation, 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, integer const useTotalMassEquation, 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.