Memosa-FVM  0.2
FlowModelPressureBC.h File Reference
This graph shows which files directly or indirectly include this file:

Go to the source code of this file.

Macros

#define _FLOWMODELPRESSUREBCS_H_
 

Functions

void fixedPressureMomentumBC (const StorageSite &faces, const Mesh &mesh, MultiFieldMatrix &matrix, const MultiField &xField, MultiField &rField)
 
fixedPressureContinuityBC (const StorageSite &faces, const Mesh &mesh, MultiFieldMatrix &matrix, const MultiField &xField, MultiField &rField)
 
void pressureBoundaryPostContinuitySolve (const StorageSite &faces, const Mesh &mesh, const T bp)
 
void updateFacePressureBoundary (const Mesh &mesh, const StorageSite &faces)
 

Macro Definition Documentation

#define _FLOWMODELPRESSUREBCS_H_

Definition at line 7 of file FlowModel_impl.h.

Function Documentation

T fixedPressureContinuityBC ( const StorageSite faces,
const Mesh mesh,
MultiFieldMatrix matrix,
const MultiField xField,
MultiField rField 
)

Definition at line 58 of file FlowModelPressureBC.h.

References dot(), Mesh::getCells(), StorageSite::getCount(), Mesh::getFaceCells(), MultiFieldMatrix::getMatrix(), and mag2().

Referenced by FlowModel< T >::Impl::linearizeContinuity().

62 {
63  const StorageSite& cells = mesh.getCells();
64  MultiField::ArrayIndex pIndex(&_flowFields.pressure,&cells);
65  MultiField::ArrayIndex vIndex(&_flowFields.velocity,&cells);
66  MultiField::ArrayIndex mfIndex(&_flowFields.massFlux,&faces);
67 
68  const VectorT3Array& faceArea =
69  dynamic_cast<const VectorT3Array&>(_geomFields.area[faces]);
70 
71  const TArray& cellVolume =
72  dynamic_cast<const TArray&>(_geomFields.volume[cells]);
73 
74  const VectorT3Array& cellCentroid =
75  dynamic_cast<const VectorT3Array&>(_geomFields.coordinate[cells]);
76 
77  const CRConnectivity& faceCells = mesh.getFaceCells(faces);
78 
79  PPMatrix& ppMatrix =
80  dynamic_cast<PPMatrix&>(matrix.getMatrix(pIndex,pIndex));
81 
82  const VVDiagArray& momAp = dynamic_cast<const VVDiagArray&>((*_momApField)[cells]);
83 
84  const VectorT3Array& V = dynamic_cast<const VectorT3Array&>(_flowFields.velocity[cells]);
85  const VectorT3Array& Vprev = dynamic_cast<const VectorT3Array&>((*_previousVelocity)[cells]);
86 
87  const TArray& p = dynamic_cast<const TArray&>(_flowFields.pressure[cells]);
88  const PGradArray& pGrad = dynamic_cast<const PGradArray&>(_flowFields.pressureGradient[cells]);
89 
90  const TArray& rho = dynamic_cast<TArray&>(_flowFields.density[cells]);
91 
92  PPAssembler& ppAssembler = ppMatrix.getPairWiseAssembler(faceCells);
93  PPDiagArray& ppDiag = ppMatrix.getDiag();
94 
95  TArray& rCell = dynamic_cast<TArray&>(rField[pIndex]);
96  TArray& massFlux = dynamic_cast<TArray&>(_flowFields.massFlux[faces]);
97 
98  FMatrix& dFluxdP = dynamic_cast<FMatrix&>(matrix.getMatrix(mfIndex,pIndex));
99 
100  const T momURF(_options["momentumURF"]);
101  const T OneMinusmomURF(T(1.0)-momURF);
102 
103  const int nFaces = faces.getCount();
104 
105  T netFlux(0.);
106 
107  for(int f=0; f<nFaces; f++)
108  {
109  const int c0 = faceCells(f,0);
110  const int c1 = faceCells(f,1);
111  const VectorT3 ds=cellCentroid[c1]-cellCentroid[c0];
112  const VectorT3& Af = faceArea[f];
113 
114  const T dpf = pGrad[c0]*ds - p[c1] + p[c0] ;
115  const T rhoF = rho[c0];
116 
117  // Q < 0
118  const T Q = rhoF*(Af[0]*Af[0] / momAp[c0][0] +
119  Af[1]*Af[1] / momAp[c0][1] +
120  Af[2]*Af[2] / momAp[c0][2])*cellVolume[c0]/(dot(Af,ds));
121 
122  const T massFluxI = rhoF*(dot(V[c0],Af) - OneMinusmomURF*dot(Vprev[c0],Af)) - Q*dpf +
123  OneMinusmomURF*massFlux[f];
124 
125  const VectorT3& Vb = V[c1];
126  const T massFluxB = rhoF*dot(Vb,Af);
127 
128  netFlux += massFluxI;
129 
130  massFlux[f] = massFluxI;
131 
132  T Vb_dpdVb(0);
133  if (massFluxB < 0)
134  Vb_dpdVb = -mag2(Vb)*rhoF;
135 
136 
137  const T denom = massFluxI-Q*Vb_dpdVb;
138 
139  if (denom != 0)
140  {
141  const T dMassFluxdp0 = -Q*massFluxI/denom;
142  dFluxdP.setCoeffL(f,dMassFluxdp0);
143  dFluxdP.setCoeffR(f,T(0.));
144 
145  const T dpbdp0 = -Q*Vb_dpdVb/denom;
146  // contribution to cell equation
147  rCell[c0] -= massFlux[f];
148  ppDiag[c0] -= dMassFluxdp0;
149  ppDiag[c1] = -1;
150  ppAssembler.getCoeff01(f) = 0.;
151  ppAssembler.getCoeff10(f) = dpbdp0;
152  }
153  else
154  {
155  // treat as fixed pressure
156  dFluxdP.setCoeffL(f,-Q);
157  dFluxdP.setCoeffR(f,T(0.));
158  ppDiag[c0] += Q;
159  ppDiag[c1] = -1;
160  rCell[c0] -= massFlux[f];
161  rCell[c1] = 0;
162  ppAssembler.getCoeff10(f) = 0.;
163  ppAssembler.getCoeff01(f) = 0.;
164  }
165  ppMatrix.setBoundary(c1);
166  }
167  return netFlux;
168 }
Matrix & getMatrix(const Index &rowIndex, const Index &colIndex)
pair< const Field *, const StorageSite * > ArrayIndex
Definition: MultiField.h:21
const StorageSite & getCells() const
Definition: Mesh.h:109
const CRConnectivity & getFaceCells(const StorageSite &site) const
Definition: Mesh.cpp:388
int getCount() const
Definition: StorageSite.h:39
T dot(const Vector< T, 3 > &a, const Vector< T, 3 > &b)
Definition: Vector.h:253
T mag2(const Vector< T, 3 > &a)
Definition: Vector.h:267
void fixedPressureMomentumBC ( const StorageSite faces,
const Mesh mesh,
MultiFieldMatrix matrix,
const MultiField xField,
MultiField rField 
)

Definition at line 11 of file FlowModelPressureBC.h.

References Mesh::getCells(), StorageSite::getCount(), Mesh::getFaceCells(), MultiFieldMatrix::getMatrix(), and mag2().

Referenced by FlowModel< T >::Impl::linearizeMomentum().

15 {
16  const StorageSite& cells = mesh.getCells();
17  MultiField::ArrayIndex vIndex(&_flowFields.velocity,&cells);
18 
19  const CRConnectivity& faceCells = mesh.getFaceCells(faces);
20 
21  const VectorT3Array& faceArea =
22  dynamic_cast<const VectorT3Array&>(_geomFields.area[faces]);
23 
24  const TArray& faceAreaMag =
25  dynamic_cast<const TArray&>(_geomFields.areaMag[faces]);
26 
27  VVMatrix& vvMatrix =
28  dynamic_cast<VVMatrix&>(matrix.getMatrix(vIndex,vIndex));
29 
30  VVDiagArray& vvDiag = vvMatrix.getDiag();
31 
32  const VectorT3Array& V = dynamic_cast<const VectorT3Array&>(_flowFields.velocity[cells]);
33 
34  const TArray& rho = dynamic_cast<TArray&>(_flowFields.density[cells]);
35 
36 
37  TArray& massFlux = dynamic_cast<TArray&>(_flowFields.massFlux[faces]);
38 
39  const T momURF(_options["momentumURF"]);
40 
41  const int nFaces = faces.getCount();
42 
43  for(int f=0; f<nFaces; f++)
44  {
45  if (massFlux[f] < 0.)
46  {
47  const int c0 = faceCells(f,0);
48  const int c1 = faceCells(f,1);
49  const VectorT3& Af = faceArea[f];
50  const T dpdV = -rho[c0]*mag2(V[c1])/momURF;
51  vvDiag[c0][0] += dpdV*Af[0]*Af[0]/faceAreaMag[f];
52  vvDiag[c0][1] += dpdV*Af[1]*Af[1]/faceAreaMag[f];
53  vvDiag[c0][2] += dpdV*Af[2]*Af[2]/faceAreaMag[f];
54  }
55  }
56 }
Matrix & getMatrix(const Index &rowIndex, const Index &colIndex)
pair< const Field *, const StorageSite * > ArrayIndex
Definition: MultiField.h:21
const StorageSite & getCells() const
Definition: Mesh.h:109
const CRConnectivity & getFaceCells(const StorageSite &site) const
Definition: Mesh.cpp:388
int getCount() const
Definition: StorageSite.h:39
T mag2(const Vector< T, 3 > &a)
Definition: Vector.h:267
void pressureBoundaryPostContinuitySolve ( const StorageSite faces,
const Mesh mesh,
const T  bp 
)

Definition at line 170 of file FlowModelPressureBC.h.

References Mesh::getCells(), StorageSite::getCount(), Mesh::getFaceCells(), and mag2().

Referenced by FlowModel< T >::Impl::postContinuitySolve().

173 {
174  const StorageSite& cells = mesh.getCells();
175 
176  const VectorT3Array& faceArea =
177  dynamic_cast<const VectorT3Array&>(_geomFields.area[faces]);
178 
179  const TArray& faceAreaMag =
180  dynamic_cast<const TArray&>(_geomFields.areaMag[faces]);
181 
182 
183  const CRConnectivity& faceCells = mesh.getFaceCells(faces);
184  VectorT3Array& V = dynamic_cast<VectorT3Array&>(_flowFields.velocity[cells]);
185  TArray& p = dynamic_cast<TArray&>(_flowFields.pressure[cells]);
186  const TArray& rho = dynamic_cast<TArray&>(_flowFields.density[cells]);
187  TArray& massFlux = dynamic_cast<TArray&>(_flowFields.massFlux[faces]);
188 
189  const int nFaces = faces.getCount();
190  for(int f=0; f<nFaces; f++)
191  {
192  const int c0 = faceCells(f,0);
193  const int c1 = faceCells(f,1);
194  const T rhoF = rho[c0];
195 
196  if (massFlux[f] > 0)
197  {
198  //const T Vn = massFlux[f]/(rhoF*faceAreaMag[f]);
199  //const VectorT3 Vt = V[c0]
200  // - dot(V[c0],faceArea[f])/(faceAreaMag[f]*faceAreaMag[f])*faceArea[f];
201  // V[c1] = Vn/faceAreaMag[f]*faceArea[f]+Vt ;
202  V[c1] = V[c0];
203  p[c1]=bp;
204  }
205  else
206  {
207  const T Vn = -massFlux[f]/(rhoF*faceAreaMag[f]);
208 
209  V[c1] = -Vn*faceArea[f]/faceAreaMag[f];
210  p[c1] = bp - 0.5*rhoF*mag2(V[c1]);
211  }
212  }
213 }
const StorageSite & getCells() const
Definition: Mesh.h:109
const CRConnectivity & getFaceCells(const StorageSite &site) const
Definition: Mesh.cpp:388
int getCount() const
Definition: StorageSite.h:39
T mag2(const Vector< T, 3 > &a)
Definition: Vector.h:267
void updateFacePressureBoundary ( const Mesh mesh,
const StorageSite faces 
)

Definition at line 218 of file FlowModelPressureBC.h.

References Mesh::getCells(), StorageSite::getCount(), and Mesh::getFaceCells().

Referenced by FlowModel< T >::Impl::postContinuitySolve().

220 {
221  const StorageSite& cells = mesh.getCells();
222 
223  const CRConnectivity& faceCells = mesh.getFaceCells(faces);
224 
225  const TArray& pCell = dynamic_cast<const TArray&>(_flowFields.pressure[cells]);
226  TArray& pFace = dynamic_cast<TArray&>(_flowFields.pressure[faces]);
227 
228  const int nFaces = faces.getCount();
229  for(int f=0; f<nFaces; f++)
230  {
231  // const int c0 = faceCells(f,0);
232  const int c1 = faceCells(f,1);
233 
234  pFace[f] = pCell[c1];
235  }
236 }
const StorageSite & getCells() const
Definition: Mesh.h:109
const CRConnectivity & getFaceCells(const StorageSite &site) const
Definition: Mesh.cpp:388
int getCount() const
Definition: StorageSite.h:39