GEOS
MultiphasePoromechanics_impl.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_MULTIPHYSICS_POROMECHANICSKERNELS_MULTIPHASEPOROMECHANICS_IMPL_HPP_
21 #define GEOS_PHYSICSSOLVERS_MULTIPHYSICS_POROMECHANICSKERNELS_MULTIPHASEPOROMECHANICS_IMPL_HPP_
22 
23 #include "constitutive/fluid/multifluid/MultiFluidBase.hpp"
24 #include "finiteElement/BilinearFormUtilities.hpp"
25 #include "finiteElement/LinearFormUtilities.hpp"
26 #include "finiteElement/elementFormulations/FiniteElementOperators.hpp"
30 
31 namespace geos
32 {
33 
34 namespace poromechanicsKernels
35 {
36 
37 template< typename SUBREGION_TYPE,
38  typename CONSTITUTIVE_TYPE,
39  typename FE_TYPE >
41 MultiphasePoromechanics( NodeManager const & nodeManager,
42  EdgeManager const & edgeManager,
43  FaceManager const & faceManager,
44  localIndex const targetRegionIndex,
45  SUBREGION_TYPE const & elementSubRegion,
46  FE_TYPE const & finiteElementSpace,
47  CONSTITUTIVE_TYPE & inputConstitutiveType,
48  arrayView1d< globalIndex const > const inputDispDofNumber,
49  globalIndex const rankOffset,
51  arrayView1d< real64 > const inputRhs,
52  real64 const inputDt,
53  real64 const (&gravityVector)[3],
54  string const inputFlowDofKey,
55  localIndex const numComponents,
56  localIndex const numPhases,
57  integer const useSimpleAccumulation,
58  integer const useTotalMassEquation,
59  integer const performStressInitialization,
60  string const fluidModelKey ):
61  Base( nodeManager,
62  edgeManager,
63  faceManager,
64  targetRegionIndex,
65  elementSubRegion,
66  finiteElementSpace,
67  inputConstitutiveType,
68  inputDispDofNumber,
69  rankOffset,
70  inputMatrix,
71  inputRhs,
72  inputDt,
73  gravityVector,
74  inputFlowDofKey,
75  fluidModelKey ),
76  m_fluidPhaseVolFrac( elementSubRegion.template getField< fields::flow::phaseVolumeFraction >() ),
77  m_fluidPhaseVolFrac_n( elementSubRegion.template getField< fields::flow::phaseVolumeFraction_n >() ),
78  m_dFluidPhaseVolFrac( elementSubRegion.template getField< fields::flow::dPhaseVolumeFraction >() ),
79  m_dGlobalCompFraction_dGlobalCompDensity( elementSubRegion.template getField< fields::flow::dGlobalCompFraction_dGlobalCompDensity >() ),
80  m_compDens( elementSubRegion.template getField< fields::flow::globalCompDensity >() ),
81  m_compDens_n( elementSubRegion.template getField< fields::flow::globalCompDensity_n >() ),
82  m_numComponents( numComponents ),
83  m_numPhases( numPhases ),
84  m_useSimpleAccumulation( useSimpleAccumulation ),
85  m_useTotalMassEquation( useTotalMassEquation ),
86  m_performStressInitialization( performStressInitialization )
87 {
88  GEOS_ERROR_IF_GT_MSG( m_numComponents, maxNumComponents,
89  "MultiphasePoromechanics solver allows at most " <<
90  maxNumComponents << " components at the moment" );
91 
92  // extract fluid constitutive data views
93  {
94  string const & fluidModelName = elementSubRegion.template getReference< string >( fluidModelKey );
95  constitutive::MultiFluidBase const & fluid =
96  elementSubRegion.template getConstitutiveModel< constitutive::MultiFluidBase >( fluidModelName );
97 
98  m_fluidPhaseDensity = fluid.phaseDensity();
99  m_fluidPhaseDensity_n = fluid.phaseDensity_n();
100  m_dFluidPhaseDensity = fluid.dPhaseDensity();
101 
102  m_fluidPhaseCompFrac = fluid.phaseCompFraction();
103  m_fluidPhaseCompFrac_n = fluid.phaseCompFraction_n();
104  m_dFluidPhaseCompFrac = fluid.dPhaseCompFraction();
105 
106  m_fluidPhaseMassDensity = fluid.phaseMassDensity();
107  m_dFluidPhaseMassDensity = fluid.dPhaseMassDensity();
108 
109  m_useMass = fluid.getMassFlag();
110  m_componentMolarWeight = fluid.componentMolarWeights();
111  }
112 }
113 
122 template< typename SUBREGION_TYPE,
123  typename CONSTITUTIVE_TYPE,
124  typename FE_TYPE >
128 setup( localIndex const k,
129  StackVariables & stack ) const
130 {
131  // initialize displacement dof and pressure dofs
132  Base::setup( k, stack );
133 
134  // setup component dofs
135  // for now, maxNumComponents > m_numComponents, so we pad localComponentDofIndices with -1
136  LvArray::tensorOps::fill< maxNumComponents >( stack.localComponentDofIndices, -1.0 );
137  for( integer flowDofIndex=0; flowDofIndex < m_numComponents; ++flowDofIndex )
138  {
139  stack.localComponentDofIndices[flowDofIndex] = stack.localPressureDofIndex + flowDofIndex + 1;
140  }
141 }
142 
143 template< typename SUBREGION_TYPE,
144  typename CONSTITUTIVE_TYPE,
145  typename FE_TYPE >
150  localIndex const q,
151  StackVariables & stack ) const
152 {
153  real64 porosity = 0.0;
154  real64 porosity_n = 0.0;
155  real64 dPorosity_dVolStrain = 0.0;
156  real64 dPorosity_dPressure = 0.0;
157  real64 dPorosity_dTemperature = 0.0;
158  real64 dSolidDensity_dPressure = 0.0;
159 
160  // Step 1: call the constitutive model to evaluate the total stress and compute porosity
161  m_constitutiveUpdate.smallStrainUpdatePoromechanics( k, q,
162  m_dt,
163  m_pressure[k],
164  m_pressure_n[k],
165  stack.deltaTemperature,
167  stack.strainIncrement,
168  stack.totalStress,
171  stack.stiffness,
172  m_performStressInitialization,
173  porosity,
174  porosity_n,
175  dPorosity_dVolStrain,
176  dPorosity_dPressure,
177  dPorosity_dTemperature,
178  dSolidDensity_dPressure );
179 
180  // Step 2: compute the body force
181  computeBodyForce( k, q,
182  porosity,
183  dPorosity_dVolStrain,
184  dPorosity_dPressure,
185  dPorosity_dTemperature,
186  dSolidDensity_dPressure,
187  stack );
188 
189  // Step 3: compute fluid mass increment
190  computeFluidIncrement( k, q,
191  porosity,
192  porosity_n,
193  dPorosity_dVolStrain,
194  dPorosity_dPressure,
195  dPorosity_dTemperature,
196  stack );
197 
198  // Step 4: compute pore volume constraint
199  computePoreVolumeConstraint( k,
200  porosity_n,
201  stack );
202 }
203 
204 template< typename SUBREGION_TYPE,
205  typename CONSTITUTIVE_TYPE,
206  typename FE_TYPE >
207 template< typename FUNC >
212  localIndex const q,
213  real64 const & porosity,
214  real64 const & dPorosity_dVolStrain,
215  real64 const & dPorosity_dPressure,
216  real64 const & dPorosity_dTemperature,
217  real64 const & dSolidDensity_dPressure,
218  StackVariables & stack,
219  FUNC && bodyForceKernelOp ) const
220 {
221  GEOS_UNUSED_VAR( dPorosity_dTemperature );
222 
223  // Step 1: compute fluid total mass density and its derivatives
224 
225  real64 totalMassDensity = 0.0;
226  real64 dTotalMassDensity_dComponents[maxNumComponents]{};
227 
228  arraySlice1d< real64 const, compflow::USD_COMP - 1 > const compDens = m_compDens[k];
229  for( integer ic = 0; ic < m_numComponents; ++ic )
230  {
231  totalMassDensity = totalMassDensity + ( m_useMass ? compDens[ic] : compDens[ic] * m_componentMolarWeight[ic] );
232  dTotalMassDensity_dComponents[ic] = m_useMass ? 1.0 : m_componentMolarWeight[ic];
233  }
234 
235  // Step 2: compute mixture density as an average between total mass density and solid density
236 
237  real64 const mixtureDensity = ( 1.0 - porosity ) * m_solidDensity( k, q ) + porosity * totalMassDensity;
238  real64 const dMixtureDens_dVolStrainIncrement = dPorosity_dVolStrain * ( -m_solidDensity( k, q ) + totalMassDensity );
239  real64 const dMixtureDens_dPressure = dPorosity_dPressure * ( -m_solidDensity( k, q ) + totalMassDensity )
240  + ( 1.0 - porosity ) * dSolidDensity_dPressure;
241  LvArray::tensorOps::scale< maxNumComponents >( dTotalMassDensity_dComponents, porosity );
242 
243  // Step 3: finally, get the body force
244 
245  LvArray::tensorOps::scaledCopy< 3 >( stack.bodyForce, m_gravityVector, mixtureDensity );
246  LvArray::tensorOps::scaledCopy< 3 >( stack.dBodyForce_dVolStrainIncrement, m_gravityVector, dMixtureDens_dVolStrainIncrement );
247  LvArray::tensorOps::scaledCopy< 3 >( stack.dBodyForce_dPressure, m_gravityVector, dMixtureDens_dPressure );
248  LvArray::tensorOps::Rij_eq_AiBj< 3, maxNumComponents >( stack.dBodyForce_dComponents, m_gravityVector, dTotalMassDensity_dComponents );
249 
250  // Step 4: customize the kernel (for instance, to add thermal derivatives)
251  bodyForceKernelOp( totalMassDensity, mixtureDensity );
252 
253 }
254 
255 template< typename SUBREGION_TYPE,
256  typename CONSTITUTIVE_TYPE,
257  typename FE_TYPE >
258 template< typename FUNC >
263  localIndex const q,
264  real64 const & porosity,
265  real64 const & porosity_n,
266  real64 const & dPorosity_dVolStrain,
267  real64 const & dPorosity_dPressure,
268  real64 const & dPorosity_dTemperature,
269  StackVariables & stack,
270  FUNC && fluidIncrementKernelOp ) const
271 {
272  LvArray::tensorOps::fill< maxNumComponents >( stack.compMassIncrement, 0.0 );
273  LvArray::tensorOps::fill< maxNumComponents >( stack.dCompMassIncrement_dVolStrainIncrement, 0.0 );
274  LvArray::tensorOps::fill< maxNumComponents >( stack.dCompMassIncrement_dPressure, 0.0 );
275  LvArray::tensorOps::fill< maxNumComponents, maxNumComponents >( stack.dCompMassIncrement_dComponents, 0.0 );
276 
277  if( m_useSimpleAccumulation )
278  {
279  arraySlice1d< real64 const, compflow::USD_COMP - 1 >
280  const compDens = m_compDens[k];
281  arraySlice1d< real64 const, compflow::USD_COMP - 1 >
282  const compDens_n = m_compDens_n[k];
283 
284  for( integer ic = 0; ic < m_numComponents; ++ic )
285  {
286  stack.compMassIncrement[ic] += porosity * compDens[ic] - porosity_n * compDens_n[ic];
287  stack.dCompMassIncrement_dPressure[ic] += dPorosity_dPressure * compDens[ic];
288  stack.dCompMassIncrement_dVolStrainIncrement[ic] += dPorosity_dVolStrain * compDens[ic];
289  stack.dCompMassIncrement_dComponents[ic][ic] += porosity;
290  }
291  }
292  else
293  {
294  using Deriv = constitutive::multifluid::DerivativeOffset;
295 
296  GEOS_UNUSED_VAR( dPorosity_dTemperature );
297 
298  // temporary work arrays and slices
299  real64 dPhaseAmount_dC[maxNumComponents]{};
300  real64 dPhaseCompFrac_dC[maxNumComponents]{};
301 
302  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 >
303  const phaseDensity = m_fluidPhaseDensity[k][q];
304  arraySlice1d< real64 const, constitutive::multifluid::USD_PHASE - 2 >
305  const phaseDensity_n = m_fluidPhaseDensity_n[k][q];
306  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_DC - 2 >
307  const dPhaseDensity = m_dFluidPhaseDensity[k][q];
308  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 >
309  const phaseCompFrac = m_fluidPhaseCompFrac[k][q];
310  arraySlice2d< real64 const, constitutive::multifluid::USD_PHASE_COMP - 2 >
311  const phaseCompFrac_n = m_fluidPhaseCompFrac_n[k][q];
312  arraySlice3d< real64 const, constitutive::multifluid::USD_PHASE_COMP_DC - 2 >
313  const dPhaseCompFrac = m_dFluidPhaseCompFrac[k][q];
314  arraySlice1d< real64 const, compflow::USD_PHASE - 1 >
315  const phaseVolFrac = m_fluidPhaseVolFrac[k];
316  arraySlice1d< real64 const, compflow::USD_PHASE - 1 >
317  const phaseVolFrac_n = m_fluidPhaseVolFrac_n[k];
318  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 >
319  const dPhaseVolFrac = m_dFluidPhaseVolFrac[k];
320  arraySlice2d< real64 const, compflow::USD_COMP_DC - 1 >
321  const dGlobalCompFrac_dGlobalCompDensity = m_dGlobalCompFraction_dGlobalCompDensity[k];
322 
323  // loop over the fluid phases
324  for( integer ip = 0; ip < m_numPhases; ++ip )
325  {
326 
327  // compute the mass of the current phase
328  real64 const phaseAmount = porosity * phaseVolFrac( ip ) * phaseDensity( ip );
329  real64 const phaseAmount_n = porosity_n * phaseVolFrac_n( ip ) * phaseDensity_n( ip );
330 
331  real64 const dPhaseAmount_dVolStrain = dPorosity_dVolStrain * phaseVolFrac( ip ) * phaseDensity( ip );
332  real64 const dPhaseAmount_dP = dPorosity_dPressure * phaseVolFrac( ip ) * phaseDensity( ip )
333  + porosity * ( dPhaseVolFrac( ip, Deriv::dP ) * phaseDensity( ip )
334  + phaseVolFrac( ip ) * dPhaseDensity( ip, Deriv::dP ) );
335 
336  applyChainRule( m_numComponents,
337  dGlobalCompFrac_dGlobalCompDensity,
338  dPhaseDensity[ip],
339  dPhaseAmount_dC,
340  Deriv::dC );
341 
342  for( integer jc = 0; jc < m_numComponents; ++jc )
343  {
344  dPhaseAmount_dC[jc] = dPhaseAmount_dC[jc] * phaseVolFrac( ip )
345  + phaseDensity( ip ) * dPhaseVolFrac( ip, Deriv::dC + jc );
346  dPhaseAmount_dC[jc] *= porosity;
347  }
348 
349  // for each phase, compute the amount of each component transported by the phase
350  for( integer ic = 0; ic < m_numComponents; ++ic )
351  {
352  stack.compMassIncrement[ic] += phaseAmount * phaseCompFrac( ip, ic )
353  - phaseAmount_n * phaseCompFrac_n( ip, ic );
354 
355  stack.dCompMassIncrement_dPressure[ic] += dPhaseAmount_dP * phaseCompFrac( ip, ic )
356  + phaseAmount * dPhaseCompFrac( ip, ic, Deriv::dP );
357  stack.dCompMassIncrement_dVolStrainIncrement[ic] += dPhaseAmount_dVolStrain * phaseCompFrac( ip, ic );
358 
359  applyChainRule( m_numComponents,
360  dGlobalCompFrac_dGlobalCompDensity,
361  dPhaseCompFrac[ip][ic],
362  dPhaseCompFrac_dC,
363  Deriv::dC );
364 
365  for( integer jc = 0; jc < m_numComponents; ++jc )
366  {
367  stack.dCompMassIncrement_dComponents[ic][jc] += dPhaseAmount_dC[jc] * phaseCompFrac( ip, ic )
368  + phaseAmount * dPhaseCompFrac_dC[jc];
369  }
370  }
371 
372  // call the lambda in the phase loop to allow the reuse of the phase amounts and their derivatives
373  // possible use: assemble the derivatives wrt temperature, and the accumulation term of the energy equation for this phase
374  fluidIncrementKernelOp( ip, phaseAmount, phaseAmount_n, dPhaseAmount_dVolStrain, dPhaseAmount_dP, dPhaseAmount_dC );
375 
376  }
377  }
378 }
379 
380 template< typename SUBREGION_TYPE,
381  typename CONSTITUTIVE_TYPE,
382  typename FE_TYPE >
387  real64 const & porosity_n,
388  StackVariables & stack ) const
389 {
390  using Deriv = constitutive::multifluid::DerivativeOffset;
391 
392  arraySlice1d< real64 const, compflow::USD_PHASE - 1 > const phaseVolFrac = m_fluidPhaseVolFrac[k];
393  arraySlice2d< real64 const, compflow::USD_PHASE_DC - 1 > const dPhaseVolFrac = m_dFluidPhaseVolFrac[k];
394 
395  stack.poreVolConstraint = 1.0;
396  stack.dPoreVolConstraint_dPressure = 0.0;
397  LvArray::tensorOps::fill< 1, maxNumComponents >( stack.dPoreVolConstraint_dComponents, 0.0 );
398 
399  for( integer ip = 0; ip < m_numPhases; ++ip )
400  {
401  stack.poreVolConstraint -= phaseVolFrac( ip );
402  stack.dPoreVolConstraint_dPressure -= dPhaseVolFrac( ip, Deriv::dP ) * porosity_n;
403 
404  for( integer jc = 0; jc < m_numComponents; ++jc )
405  {
406  stack.dPoreVolConstraint_dComponents[0][jc] -= dPhaseVolFrac( ip, Deriv::dC+jc ) * porosity_n;
407  }
408  }
409  stack.poreVolConstraint *= porosity_n;
410 }
411 
412 template< typename SUBREGION_TYPE,
413  typename CONSTITUTIVE_TYPE,
414  typename FE_TYPE >
418 assembleMomentumBalanceTerms( real64 const ( &N )[numNodesPerElem],
419  real64 const ( &dNdX )[numNodesPerElem][3],
420  real64 const & detJxW,
421  StackVariables & stack ) const
422 {
423  using namespace PDEUtilities;
424 
425  constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
426  constexpr FunctionSpace displacementTestSpace = displacementTrialSpace;
427  constexpr FunctionSpace pressureTrialSpace = FunctionSpace::P0;
428 
429  // Step 1: compute local linear momentum balance residual
430 
431  LinearFormUtilities::compute< displacementTestSpace,
432  DifferentialOperator::SymmetricGradient >
433  (
434  stack.localResidualMomentum,
435  dNdX,
436  stack.totalStress,
437  -detJxW );
438 
439  LinearFormUtilities::compute< displacementTestSpace,
440  DifferentialOperator::Identity >
441  (
442  stack.localResidualMomentum,
443  N,
444  stack.bodyForce,
445  detJxW );
446 
447  // Step 2: compute local linear momentum balance residual derivatives with respect to displacement
448  BilinearFormUtilities::compute< displacementTestSpace,
449  displacementTrialSpace,
450  DifferentialOperator::SymmetricGradient,
451  DifferentialOperator::SymmetricGradient >
452  (
454  dNdX,
455  stack.stiffness, // fourth-order tensor handled via DiscretizationOps
456  dNdX,
457  -detJxW );
458 
459  BilinearFormUtilities::compute< displacementTestSpace,
460  displacementTrialSpace,
461  DifferentialOperator::Identity,
462  DifferentialOperator::Divergence >
463  (
465  N,
467  dNdX,
468  detJxW );
469 
470  // Step 3: compute local linear momentum balance residual derivatives with respect to pressure
471  BilinearFormUtilities::compute< displacementTestSpace,
472  pressureTrialSpace,
473  DifferentialOperator::SymmetricGradient,
474  DifferentialOperator::Identity >
475  (
477  dNdX,
479  1.0,
480  -detJxW );
481 
482  BilinearFormUtilities::compute< displacementTestSpace,
483  pressureTrialSpace,
484  DifferentialOperator::Identity,
485  DifferentialOperator::Identity >
486  (
488  N,
489  stack.dBodyForce_dPressure,
490  1.0,
491  detJxW );
492 
493  // Step 4: compute local linear momentum balance residual derivatives with respect to components
494  BilinearFormUtilities::compute< displacementTestSpace,
495  FunctionSpace::P0,
496  DifferentialOperator::Identity,
497  DifferentialOperator::Identity >
498  (
500  N,
502  1.0,
503  detJxW );
504 }
505 
506 template< typename SUBREGION_TYPE,
507  typename CONSTITUTIVE_TYPE,
508  typename FE_TYPE >
512 assembleElementBasedFlowTerms( real64 const ( &dNdX )[numNodesPerElem][3],
513  real64 const & detJxW,
514  StackVariables & stack ) const
515 {
516  using namespace PDEUtilities;
517 
518  constexpr FunctionSpace displacementTrialSpace = FE_TYPE::template getFunctionSpace< numDofPerTrialSupportPoint >();
519  constexpr FunctionSpace displacementTestSpace = displacementTrialSpace;
520 
521  // Step 1: mass balance equations
522 
523  // compute local component mass balance residual
524  LinearFormUtilities::compute< FunctionSpace::P0,
525  DifferentialOperator::Identity >
526  (
527  stack.localResidualMass,
528  1.0,
529  stack.compMassIncrement,
530  detJxW );
531 
532  // compute local mass balance residual derivatives with respect to displacement
533  BilinearFormUtilities::compute< FunctionSpace::P0,
534  displacementTestSpace,
535  DifferentialOperator::Identity,
536  DifferentialOperator::Divergence >
537  (
539  1.0,
541  dNdX,
542  detJxW );
543 
544  // compute local mass balance residual derivatives with respect to pressure
545  BilinearFormUtilities::compute< FunctionSpace::P0,
546  FunctionSpace::P0,
547  DifferentialOperator::Identity,
548  DifferentialOperator::Identity >
549  (
551  1.0,
553  1.0,
554  detJxW );
555 
556  // compute local mass balance residual derivatives with respect to components
557  BilinearFormUtilities::compute< FunctionSpace::P0,
558  FunctionSpace::P0,
559  DifferentialOperator::Identity,
560  DifferentialOperator::Identity >
561  (
563  1.0,
565  1.0,
566  detJxW );
567 
568 
569  // Step 2: pore volume constraint equation
570 
571  // compute local pore volume contraint residual
572  LinearFormUtilities::compute< FunctionSpace::P0,
573  DifferentialOperator::Identity >
574  (
576  1.0,
577  stack.poreVolConstraint,
578  detJxW );
579 
580  // compute local pore volume contraint residual derivatives with respect to pressure
581  BilinearFormUtilities::compute< FunctionSpace::P0,
582  FunctionSpace::P0,
583  DifferentialOperator::Identity,
584  DifferentialOperator::Identity >
585  (
587  1.0,
589  1.0,
590  detJxW );
591 
592  // compute local pore volume contraint residual derivatives with respect to components
593  BilinearFormUtilities::compute< FunctionSpace::P0,
594  FunctionSpace::P0,
595  DifferentialOperator::Identity,
596  DifferentialOperator::Identity >
597  (
599  1.0,
601  1.0,
602  detJxW );
603 }
604 
605 
606 template< typename SUBREGION_TYPE,
607  typename CONSTITUTIVE_TYPE,
608  typename FE_TYPE >
613  localIndex const q,
614  StackVariables & stack ) const
615 {
616  // Step 1: compute displacement finite element basis functions (N), basis function derivatives (dNdX), and
617  // determinant of the Jacobian transformation matrix times the quadrature weight (detJxW)
618  real64 N[numNodesPerElem]{};
619  real64 dNdX[numNodesPerElem][3]{};
620  FE_TYPE::calcN( q, stack.feStack, N );
621  real64 const detJxW = FE_TYPE::calcGradN( q, stack.xLocal,
622  stack.feStack, dNdX );
623 
624  // Step 2: compute strain increment
625  LvArray::tensorOps::fill< 6 >( stack.strainIncrement, 0.0 );
626  finiteElement::feOps::symmetricGradient( dNdX, stack.uhat_local, stack.strainIncrement );
627 
628  // Step 3: compute 1) the total stress, 2) the body force terms, and 3) the fluidMassIncrement
629  // using quantities returned by the PorousSolid constitutive model.
630  // This function also computes the derivatives of these three quantities wrt primary variables
631  smallStrainUpdate( k, q, stack );
632 
633  // Step 4: use the total stress and the body force to increment the local momentum balance residual
634  // This function also fills the local Jacobian rows corresponding to the momentum balance.
635  assembleMomentumBalanceTerms( N, dNdX, detJxW, stack );
636 
637  // Step 5: use the fluid mass increment to increment the local mass balance residual
638  // This function also fills the local Jacobian rows corresponding to the mass balance.
639  assembleElementBasedFlowTerms( dNdX, detJxW, stack );
640 }
641 
645 template< typename SUBREGION_TYPE,
646  typename CONSTITUTIVE_TYPE,
647  typename FE_TYPE >
651 complete( localIndex const k,
652  StackVariables & stack ) const
653 {
654  using namespace compositionalMultiphaseUtilities;
655 
656  GEOS_UNUSED_VAR( k );
657 
658  real64 maxForce = 0;
659  localIndex const numSupportPoints =
660  m_finiteElementSpace.getNumSupportPoints( stack.feStack );
661  integer numDisplacementDofs = numSupportPoints * numDofPerTestSupportPoint;
662  constexpr integer maxNumDisplacementDofs = FE_TYPE::maxSupportPoints * numDofPerTestSupportPoint;
663 
664  if( m_useTotalMassEquation > 0 )
665  {
666  // Apply equation/variable change transformation(s)
667  real64 work[maxNumDisplacementDofs > ( maxNumComponents + 1 ) ? maxNumDisplacementDofs : maxNumComponents + 1];
668  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( m_numComponents, numDisplacementDofs, stack.dLocalResidualMass_dDisplacement, work );
669  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( m_numComponents, 1, stack.dLocalResidualMass_dPressure, work );
670  shiftRowsAheadByOneAndReplaceFirstRowWithColumnSum( m_numComponents, m_numComponents, stack.dLocalResidualMass_dComponents, work );
671  shiftElementsAheadByOneAndReplaceFirstElementWithSum( m_numComponents, stack.localResidualMass );
672  }
673 
674  for( int localNode = 0; localNode < numSupportPoints; ++localNode )
675  {
676  for( int dim = 0; dim < numDofPerTestSupportPoint; ++dim )
677  {
678  localIndex const dof = LvArray::integerConversion< localIndex >( stack.localRowDofIndex[numDofPerTestSupportPoint * localNode + dim] - m_dofRankOffset );
679 
680  // we need this check to filter out ghost nodes in the assembly
681  if( dof < 0 || dof >= m_matrix.numRows() )
682  {
683  continue;
684  }
685  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
686  stack.localRowDofIndex,
687  stack.dLocalResidualMomentum_dDisplacement[numDofPerTestSupportPoint * localNode + dim],
688  numDisplacementDofs );
689 
690  RAJA::atomicAdd< parallelDeviceAtomic >( &m_rhs[dof], stack.localResidualMomentum[numDofPerTestSupportPoint * localNode + dim] );
691  maxForce = fmax( maxForce, fabs( stack.localResidualMomentum[numDofPerTestSupportPoint * localNode + dim] ) );
692 
693  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
694  &stack.localPressureDofIndex,
695  stack.dLocalResidualMomentum_dPressure[numDofPerTestSupportPoint * localNode + dim],
696  1 );
697 
698  m_matrix.template addToRowBinarySearchUnsorted< parallelDeviceAtomic >( dof,
700  stack.dLocalResidualMomentum_dComponents[numDofPerTestSupportPoint * localNode + dim],
701  m_numComponents );
702  }
703  }
704 
705  localIndex const dof = LvArray::integerConversion< localIndex >( stack.localPressureDofIndex - m_dofRankOffset );
706 
707  // we need this check to filter out ghost cells in the assembly
708  if( 0 <= dof && dof < m_matrix.numRows() )
709  {
710  for( localIndex i = 0; i < m_numComponents; ++i )
711  {
712  m_matrix.template addToRowBinarySearchUnsorted< serialAtomic >( dof + i,
713  stack.localRowDofIndex,
715  numDisplacementDofs );
716  m_matrix.template addToRow< serialAtomic >( dof + i,
717  &stack.localPressureDofIndex,
719  1 );
720  m_matrix.template addToRow< serialAtomic >( dof + i,
723  m_numComponents );
724  RAJA::atomicAdd< serialAtomic >( &m_rhs[dof+i], stack.localResidualMass[i] );
725  }
726 
727  m_matrix.template addToRow< serialAtomic >( dof + m_numComponents,
728  &stack.localPressureDofIndex,
730  1 );
731 
732  m_matrix.template addToRow< serialAtomic >( dof + m_numComponents,
735  m_numComponents );
736 
737  RAJA::atomicAdd< serialAtomic >( &m_rhs[dof+m_numComponents], stack.localResidualPoreVolConstraint[0] );
738  }
739  return maxForce;
740 }
741 
742 template< typename SUBREGION_TYPE,
743  typename CONSTITUTIVE_TYPE,
744  typename FE_TYPE >
745 template< typename POLICY,
746  typename KERNEL_TYPE >
748 kernelLaunch( localIndex const numElems,
749  KERNEL_TYPE const & kernelComponent )
750 {
752 
753  // Define a RAJA reduction variable to get the maximum residual contribution.
754  RAJA::ReduceMax< ReducePolicy< POLICY >, real64 > maxResidual( 0 );
755 
756  forAll< POLICY >( numElems,
757  [=] GEOS_HOST_DEVICE ( localIndex const k )
758  {
759  typename KERNEL_TYPE::StackVariables stack;
760 
761  kernelComponent.setup( k, stack );
762  for( integer q=0; q<KERNEL_TYPE::numQuadraturePointsPerElem; ++q )
763  {
764  kernelComponent.quadraturePointKernel( k, q, stack );
765  }
766  maxResidual.max( kernelComponent.complete( k, stack ) );
767  } );
768  return maxResidual.get();
769 }
770 
771 
772 
773 } // namespace poromechanicsKernels
774 
775 } // namespace geos
776 
777 #endif // GEOS_PHYSICSSOLVERS_MULTIPHYSICS_MULTIPHASEPOROMECHANICS_IMPL_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_FORCE_INLINE
Marks a function or lambda for inlining.
Definition: GeosxMacros.hpp:51
#define GEOS_ERROR_IF_GT_MSG(lhs, rhs,...)
Raise a hard error if one value compares greater than the other.
Definition: Logger.hpp:423
#define GEOS_MARK_FUNCTION
Mark function with both Caliper and NVTX if enabled.
This class provides an interface to ObjectManagerBase in order to manage edge data.
Definition: EdgeManager.hpp:43
The FaceManager class provides an interface to ObjectManagerBase in order to manage face data.
Definition: FaceManager.hpp:44
The NodeManager class provides an interface to ObjectManagerBase in order to manage node data.
Definition: NodeManager.hpp:46
Implements kernels for solving quasi-static multiphase poromechanics.
arrayView4d< real64 const, constitutive::multifluid::USD_PHASE_COMP > m_fluidPhaseCompFrac
Views on phase component fractions and derivatives.
GEOS_HOST_DEVICE void assembleElementBasedFlowTerms(real64 const (&dNdX)[numNodesPerElem][3], real64 const &detJxW, StackVariables &stack) const
Assemble the local mass balance residual and derivatives using fluid mass/energy increment.
arrayView3d< real64 const, constitutive::multifluid::USD_PHASE > m_fluidPhaseDensity
Views on phase densities and derivatives.
GEOS_HOST_DEVICE void computePoreVolumeConstraint(localIndex const k, real64 const &porosity_n, StackVariables &stack) const
Helper function to compute the pore-volume constraint and its derivatives wrt primary variables.
GEOS_HOST_DEVICE void computeBodyForce(localIndex const k, localIndex const q, real64 const &porosity, real64 const &dPorosity_dVolStrain, real64 const &dPorosity_dPressure, real64 const &dPorosity_dTemperature, real64 const &dSolidDensity_dPressure, StackVariables &stack, FUNC &&bodyForceKernelOp=NoOpFunc{}) const
Helper function to compute the body force term (\rho g) and its derivatives wrt primary variables.
arrayView3d< real64 const, constitutive::multifluid::USD_PHASE > m_fluidPhaseMassDensity
Views on phase mass densities.
static real64 kernelLaunch(localIndex const numElems, KERNEL_TYPE const &kernelComponent)
Kernel Launcher.
GEOS_HOST_DEVICE void smallStrainUpdate(localIndex const k, localIndex const q, StackVariables &stack) const
Helper function to compute 1) the total stress, 2) the body force term, and 3) the fluidMassIncrement...
GEOS_HOST_DEVICE void setup(localIndex const k, StackVariables &stack) const
Copy global values from primary field to a local stack array.
MultiphasePoromechanics(NodeManager const &nodeManager, EdgeManager const &edgeManager, FaceManager const &faceManager, localIndex const targetRegionIndex, SUBREGION_TYPE const &elementSubRegion, FE_TYPE const &finiteElementSpace, CONSTITUTIVE_TYPE &inputConstitutiveType, arrayView1d< globalIndex const > const inputDispDofNumber, globalIndex const rankOffset, CRSMatrixView< real64, globalIndex const > const inputMatrix, arrayView1d< real64 > const inputRhs, real64 const inputDt, real64 const (&gravityVector)[3], string const inputFlowDofKey, localIndex const numComponents, localIndex const numPhases, integer const useSimpleAccumulation, integer const useTotalMassEquation, integer const performStressInitialization, string const fluidModelKey)
Constructor.
GEOS_HOST_DEVICE void assembleMomentumBalanceTerms(real64 const (&N)[numNodesPerElem], real64 const (&dNdX)[numNodesPerElem][3], real64 const &detJxW, StackVariables &stack) const
Assemble the local linear momentum balance residual and derivatives using total stress and body force...
GEOS_HOST_DEVICE real64 complete(localIndex const k, StackVariables &stack) const
Performs the complete phase for the kernel.
GEOS_HOST_DEVICE void computeFluidIncrement(localIndex const k, localIndex const q, real64 const &porosity, real64 const &porosity_n, real64 const &dPorosity_dVolStrain, real64 const &dPorosity_dPressure, real64 const &dPorosity_dTemperature, StackVariables &stack, FUNC &&fluidIncrementKernelOp=NoOpFunc{}) const
Helper function to compute the fluid mass increment and its derivatives wrt primary variables.
Defines the kernel structure for solving quasi-static poromechanics.
ArrayView< T, 1 > arrayView1d
Alias for 1D array view.
Definition: DataTypes.hpp:179
GEOS_GLOBALINDEX_TYPE globalIndex
Global index type (for indexing objects across MPI partitions).
Definition: DataTypes.hpp:87
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
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
int integer
Signed integer type.
Definition: DataTypes.hpp:81
Kernel variables (dof numbers, jacobian and residual) located on the stack.
real64 dBodyForce_dComponents[3][maxNumComponents]
Derivative of body force wrt comp density.
real64 dLocalResidualMass_dPressure[maxNumComponents][1]
Derivative of mass balance residual wrt pressure.
real64 dCompMassIncrement_dVolStrainIncrement[maxNumComponents]
Derivative of mass accumulation wrt volumetric strain increment.
real64 dCompMassIncrement_dPressure[maxNumComponents]
Derivative of mass accumulation wrt pressure.
real64 dPoreVolConstraint_dComponents[1][maxNumComponents]
Derivative of pore volume constraint wrt comp density.
real64 dLocalResidualMass_dDisplacement[maxNumComponents][Base::StackVariables::numDispDofPerElem]
Derivative of mass balance residual wrt volumetric strain increment.
real64 dCompMassIncrement_dComponents[maxNumComponents][maxNumComponents]
Derivative of mass accumulation wrt comp density.
real64 dLocalResidualPoreVolConstraint_dPressure[1][1]
Derivative of pore volume constraint residual wrt pressure.
globalIndex localComponentDofIndices[maxNumComponents]
C-array storage for the element local row degrees of freedom.
real64 dLocalResidualMomentum_dComponents[Base::StackVariables::numDispDofPerElem][maxNumComponents]
Derivative of momemtum balance residual wrt comp density.
real64 localResidualPoreVolConstraint[1]
Pore volume constraint residual (sum of saturations = 1)
real64 dLocalResidualPoreVolConstraint_dComponents[1][maxNumComponents]
Derivative of pore volume constraint residual wrt components.
real64 dPoreVolConstraint_dPressure
Derivative of pore volume constraint wrt pressure.
real64 dLocalResidualMass_dComponents[maxNumComponents][maxNumComponents]
Derivative of mass balance residual wrt components.
real64(& dLocalResidualMomentum_dDisplacement)[numDispDofPerElem][numDispDofPerElem]
Derivative of linear momentum balance residual wrt displacement.
globalIndex localPressureDofIndex
C-array storage for the element local row degrees of freedom.
real64 dBodyForce_dPressure[3]
Derivative of body force wrt pressure.
real64 dTotalStress_dTemperature[6]
Derivative of total stress wrt temperature.
CONSTITUTIVE_TYPE::KernelWrapper::DiscretizationOps stiffness
Derivative of total stress wrt displacement.
real64 dLocalResidualMomentum_dPressure[numDispDofPerElem][1]
Derivative of linear momentum balance residual wrt pressure.
real64 dTotalStress_dPressure[6]
Derivative of total stress wrt pressure.
real64 dBodyForce_dVolStrainIncrement[3]
Derivative of body force wrt volumetric strain increment.
real64(& localResidualMomentum)[numDispDofPerElem]
Linear momentum balance residual.
real64 deltaTemperature
Delta temperature from reference state.
real64 deltaTemperatureFromLastStep
Delta temperature since last time step.