Rocstar  1.0
Rocstar multiphysics simulation application
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
FluRegion.H
Go to the documentation of this file.
1 /* *******************************************************************
2  * Rocstar Simulation Suite *
3  * Copyright@2015, Illinois Rocstar LLC. All rights reserved. *
4  * *
5  * Illinois Rocstar LLC *
6  * Champaign, IL *
7  * www.illinoisrocstar.com *
8  * sales@illinoisrocstar.com *
9  * *
10  * License: See LICENSE file in top level of distribution package or *
11  * http://opensource.org/licenses/NCSA *
12  *********************************************************************/
13 /* *******************************************************************
14  * THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, *
15  * EXPRESS OR IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES *
16  * OF MERCHANTABILITY, FITNESS FOR A PARTICULAR PURPOSE AND *
17  * NONINFRINGEMENT. IN NO EVENT SHALL THE CONTRIBUTORS OR *
18  * COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER *
19  * LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, *
20  * ARISING FROM, OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE *
21  * USE OR OTHER DEALINGS WITH THE SOFTWARE. *
22  *********************************************************************/
23 #ifndef _FLUREGION_H_
24 #define _FLUREGION_H_
25 
27 public:
28  double _current_time;
29  double _resid;
30  bool _unsteady;
31  vector<double> _rhof; // Density
32  vector<double> _rhovf; // Momentum density
33  vector<double> _rhoEf; // Energy density
34  vector<double> _pf; // Pressure
35  vector<double> _Tf; // Temperature
36  vector<double> _af; // Speed of sound
37  vector<double> _disp; // nodal displacements
38  vector<double> _gsp; // Speed of each face
39  FluVolumeSoln(unsigned int ncells=0,unsigned int nfaces=0,
40  unsigned int nnodes=0)
41  {
42  Resize(ncells,nfaces,nnodes);
43  };
44  void
45  Resize(unsigned int ncells=0,unsigned int nfaces=0,unsigned int nnodes=0)
46  {
47  _rhof.resize(ncells,0.0);
48  _rhovf.resize(3*ncells,0.0);
49  _rhoEf.resize(ncells,0.0);
50  _pf.resize(ncells,0.0);
51  _Tf.resize(ncells,0.0);
52  _af.resize(ncells,0.0);
53  _disp.resize(nnodes*3,0.0);
54  _gsp.resize(nfaces,0.0);
55  _unsteady = false;
56  _resid = 0.0;
57  _current_time = 0.0;
58  };
59 };
60 
61 class FluSurfSoln {
62 public:
63  vector<double> _du_alp; // Vertex displacement
64  vector<double> _rhofvf_alp; // Momentum flux
65  vector<double> _nf_alp; // face normal
66  vector<double> _rhof_alp; // density
67  vector<double> _pf; // pressure
68  vector<double> _qc; // convective heat flux
69  vector<double> _qr; // radiative heat flux
70  vector<double> _tf; // traction
71  vector<double> _Tb_alp; // Boundary temperature
72  vector<double> _mdot_alp; // Injection mass flux
73  vector<double> _Tflm_alp; // Flame temperature
74  vector<double> _gsp; // Grid Speeds
75  vector<int> _bflag; // Burning flag
76 
77  FluSurfSoln(unsigned int nfaces=0,unsigned int nvert=0)
78  {
79  Resize(nfaces,nvert);
80  };
81  void Resize(unsigned int nfaces=0,unsigned int nvert=0)
82  {
83  _gsp.resize(nfaces,0);
84  _du_alp.resize(nvert*3,0.0);
85  _rhofvf_alp.resize(3*nfaces,0.0);
86  _nf_alp.resize(3*nfaces,0.0);
87  _rhof_alp.resize(nfaces,0.0);
88  _pf.resize(nfaces,0.0);
89  _qc.resize(nfaces,0.0);
90  _qr.resize(nfaces,0.0);
91  _tf.resize(3*nfaces,0.0);
92  _Tb_alp.resize(nfaces,0.0);
93  _mdot_alp.resize(nfaces,0.0);
94  _Tflm_alp.resize(nfaces,0.0);
95  _bflag.resize(nfaces,0);
96  };
97 };
98 
99 class FluPatch : public GEM_DomainBoundary {
100 public:
101  // unsigned int _patchNo; = _id
102  unsigned int _ntri;
103  unsigned int _nquad;
105  bool InitSurfaceSoln(const string &prefix="",bool=false);
106  //#ifdef _ROCSTAR_X_
107  unsigned int _bcflag; // 0:nonburning, 1:burning, 2:noninteracting
109  bool RegisterSoln(const string &,bool);
110  //#endif
111 };
112 
114 public:
115  unsigned int _rbid;
116 };
117 
118 
119 class FluRegion : public GEM_Partition {
120 public:
121  string _casename;
122  string _casepath;
123  //#ifdef _ROCSTAR_X_
124  // string _rocstar_outpath;
125  //#endif
126  unsigned int _nnodes;
127  unsigned int _ntet;
128  unsigned int _nhex;
129  unsigned int _npyr;
130  unsigned int _npris;
131  unsigned int _npatches_total;
132  vector<FluPatch> _patches;
133  vector<FluBorder> _borders;
135  vector<pair<unsigned int,unsigned int> > _cellmap;
136  void report();
137  unsigned int nelem();
138  bool BuildCellMapping();
139  bool Cell2Elem(unsigned int, unsigned
140  int&, unsigned int&);
141  bool BuildPatchMapping(map<unsigned int,unsigned int> &,const string &p="");
142  bool WriteFluNative(const string &p="");
143  bool WriteFluCMP(const string &p="");
144  bool ReadFluCOM(const string &p="");
145  bool WriteFluCOM(const string &p="");
146  bool ReadFluDIM(const string &p="",double=0.0,bool=true);
147  bool WriteFluDIM(const string &p="",double=0.0,bool=true);
148  bool ReadFluGridASCII(const std::string&p="", double=0.0,bool=true);
149  bool WriteFluGridASCII(const string &pre="", double=0.0,bool=true);
150  bool ReadFluSolnASCII(const std::string&p="",unsigned int=0,
151  double=0.0,bool=true);
152  bool WriteFluSolnASCII(const std::string&p="",unsigned int=0,
153  double=0.0,bool=true);
154  bool WriteGSPASCII(const std::string&p="",unsigned int=0,
155  double=0.0,bool=true);
156  bool ReadGSPASCII(const std::string&p="",unsigned int=0,
157  double=0.0,bool=true);
158  bool ReadRegionASCII(const string &p="",unsigned int=0,unsigned int=0,
159  double=0.0,bool=true);
160  bool InitSurfaceSoln(const string &p="");
161  bool InitVolumeSoln(const string &p="");
162  bool CreateRegionMapFile(const string &p="",unsigned int=0,unsigned int=0);
163  void PopRemBordIndMPI();
164  bool PopRemBordIndFILE(const string &p="",double=0.0,bool=true);
165  bool ComRemeshInitData(const string &wname,double *cell_data,int nval_cells,
166  double *node_data,int nval_nodes);
167  bool ReadControlFile();
168  // bool validate_comm_list(int ncsend,int ncrecv,int *csend,int *crecv);
169  // bool PopulateFluVolumeSoln(double *data);
171  : GEM_Partition()
172  {
173  _casename = "";
174  _nnodes = 0;
175  _ntet = 0;
176  _npyr = 0;
177  _npris = 0;
178  _nhex = 0;
179  _npatches_total = 0;
180  _patches.resize(0);
181  _borders.resize(0);
182  _cell_ordering[0] = 1;
183  _cell_ordering[1] = 4;
184  _cell_ordering[2] = 3;
185  _cell_ordering[3] = 2;
186  };
188  : GEM_Partition(gp)
189  {
190  _casename = "";
191  _nnodes = gp._nc.size()/3;
192  _ntet = gp._tetconn.size()/4;
193  _npyr = gp._pyrconn.size()/5;
194  _npris = gp._prisconn.size()/6;
195  _nhex = gp._hexconn.size()/8;
196  _cell_ordering[0] = 1;
197  _cell_ordering[1] = 4;
198  _cell_ordering[2] = 3;
199  _cell_ordering[3] = 2;
200  };
201  template<typename DB>
202  bool PopulateFluPatches(const string &pre,vector<DB> &db)
203  {
204  map<unsigned int,unsigned int> patch_mapping;
205  BuildPatchMapping(patch_mapping,pre);
206  MapDomainBoundaries(db,patch_mapping,_patches);
207  vector<FluPatch>::iterator fpi = _patches.begin();
208  while(fpi != _patches.end()){
209  fpi->_ntri = fpi->_triconn.size()/3;
210  fpi->_nquad = fpi->_quadconn.size()/4;
211  fpi++;
212  }
213  return(true);
214  };
215 
216  template<typename PB>
217  bool PopulateFluBorders(vector<PB> &pb)
218  {
219  unsigned int nborders = pb.size();
220  _borders.resize(nborders);
221  unsigned int border = 0;
222  while(border < nborders){
223  _borders[border]._rpart = pb[border]._rpart;
224  _borders[border]._sendcells = pb[border]._sendcells;
225  _borders[border]._recvcells = pb[border]._recvcells;
226  _borders[border]._sharenodes = pb[border]._sharenodes;
227  _borders[border]._sendnodes = pb[border]._sendnodes;
228  _borders[border]._recvnodes = pb[border]._recvnodes;
229  _borders[border]._rbid = pb[border]._rbid;
230  _borders[border]._out = pb[border]._out;
231  border++;
232  }
233  return(true);
234  };
235  bool
236  validate_comm_list(int ncsend,int ncrecv,int *csend,int *crecv)
237  {
238  int index = 0;
239  int nreal_cell = _tetconn.size()/4 + _prisconn.size()/6 +
240  _pyrconn.size()/5 + _hexconn.size()/8 - (_ngtet + _ngpris +
241  _ngpyr + _nghex);
242  bool rval = true;
243  while(index < ncsend){
244  int ind = index++;
245  if(!(csend[ind] <= nreal_cell)){
246  if(_out)
247  *_out << "SEND CELL " << index << " is a ghost cell!!" << endl;
248  rval = false;
249  }
250  if(!(csend[ind] > 0)){
251  if(_out)
252  *_out << "SEND CELL " << index << " is zero or negative!" << endl;
253  rval = false;
254  }
255  }
256  index = 0;
257  list<int> recvcell_list;
258  while(index < ncrecv) {
259  int ind = index++;
260  if(!(crecv[ind] > nreal_cell)){
261  if(_out)
262  *_out << "RECV CELL " << index << " is a real cell!!" << endl;
263  rval = false;
264  }
265  if(!(crecv[ind] > 0)){
266  if(_out)
267  *_out << "RECV CELL " << index << " is zero or negative!" << endl;
268  rval = false;
269  }
270  bool duped = false;
271  list<int>::iterator rci = recvcell_list.begin();
272  while(rci != recvcell_list.end() && !duped){
273  if(crecv[ind] == *rci++){
274  if(_out)
275  *_out << "RECV CELL " << index
276  << " is duplicated in the receive list!"
277  << endl;
278  duped = true;
279  }
280  }
281  if(!duped)
282  recvcell_list.push_back(crecv[ind]);
283  }
284  return(rval);
285  };
286  void AddFluBorder(int rpid,int rpin,int nnshare, int nnsend,
287  int nnrecv,int ncsend,int ncrecv,
288  int *nshared,int *nsend,int *nrecv,
289  int *csend,int *crecv){
290  assert(rpid > 0 && rpin >= 0);
291  if(_debug && _out)
292  *_out << "FluRegion(" << _id << ")::AddFluBorder: Adding Border with"
293  << " region " << rpid << "." << endl;
294  FluBorder new_border;
295  new_border._rbid = rpin;
296  new_border._out = _out;
297  if(!validate_comm_list(ncsend,ncrecv,csend,crecv)){
298  if(_out)
299  *_out << "FluRegion(" << _id << ")::AddFluBorder: Validation of "
300  << "communication arrays failed, aborting." << endl;
301  exit(-1);
302  }
303  new_border.populate(rpid,nnshare,nnsend,nnrecv,ncsend,ncrecv,
304  nshared,nsend,nrecv,csend,crecv);
305  _borders.push_back(new_border);
306  };
307  void AddFluPatch(int patch_id,int ntri, int ngtri, int *tris,
308  int nquad,int ngquad, int *quads)
309  {
310  assert(ntri >= ngtri && nquad >= ngquad);
311  if(_debug && _out)
312  *_out << "FluRegion(" << _id << ")::AddFluPatch: Adding Patch with"
313  << " id " << patch_id << "." << endl;
314  FluPatch new_patch;
315  new_patch._ntri = ntri;
316  new_patch._nquad = nquad;
317  new_patch._id = patch_id;
318  new_patch._ngtri = ngtri;
319  new_patch._ngquad = ngquad;
320  int indy = 0;
321  new_patch._triconn.resize(3*ntri);
322  new_patch._quadconn.resize(4*nquad);
323  new_patch._out = _out;
324  while(indy < 3*ntri){
325  assert(tris[indy] != 0);
326  new_patch._triconn[indy] = tris[indy];
327  assert(new_patch._triconn[indy] != 0);
328  indy++;
329  }
330  indy = 0;
331  while(indy < 4*nquad){
332  assert(quads[indy] != 0);
333  new_patch._quadconn[indy] = quads[indy];
334  indy++;
335  }
336  _patches.push_back(new_patch);
337  };
338  //#ifdef _ROCSTAR_X_
339  bool RegisterFluSurfaceMesh();
340  bool RegisterSurfaceSoln(bool);
341  bool RegisterVolumeSoln(bool);
342  //#endif
343 };
344 
345 istream &SkipLines(istream &,unsigned int = 1);
346 
347 #endif
348 
349 
350 
351 
352 
353 
unsigned int _npyr
Definition: FluRegion.H:129
bool PopRemBordIndFILE(const string &p="", double=0.0, bool=true)
Definition: FluRegion.C:581
unsigned int _rbid
Definition: FluRegion.H:115
bool WriteFluCOM(const string &p="")
Definition: FluRegion.C:130
bool ReadFluSolnASCII(const std::string &p="", unsigned int=0, double=0.0, bool=true)
Definition: FluRegion.C:1229
void PopRemBordIndMPI()
Definition: FluRegion.C:1615
vector< double > _pf
Definition: FluRegion.H:34
vector< double > _du_alp
Definition: FluRegion.H:63
bool PopulateFluBorders(vector< PB > &pb)
Definition: FluRegion.H:217
vector< double > _rhoEf
Definition: FluRegion.H:33
bool ComRemeshInitData(const string &wname, double *cell_data, int nval_cells, double *node_data, int nval_nodes)
Definition: FluRegion.C:1653
unsigned int _bcflag
Definition: FluRegion.H:107
void MapDomainBoundaries(std::map< unsigned int, unsigned int > &bcmap)
Definition: GEM.C:1385
bool WriteFluCMP(const string &p="")
Definition: FluRegion.C:697
bool InitSurfaceSoln(const string &prefix="", bool=false)
Definition: FluRegion.C:882
bool WriteFluSolnASCII(const std::string &p="", unsigned int=0, double=0.0, bool=true)
Definition: FluRegion.C:1298
bool CreateRegionMapFile(const string &p="", unsigned int=0, unsigned int=0)
Definition: FluRegion.C:328
vector< double > _Tf
Definition: FluRegion.H:35
std::vector< unsigned int > _tetconn
Definition: GEM.H:295
bool ReadRegionASCII(const string &p="", unsigned int=0, unsigned int=0, double=0.0, bool=true)
Definition: FluRegion.C:1577
bool ReadControlFile()
Definition: FluRegion.C:101
void AddFluBorder(int rpid, int rpin, int nnshare, int nnsend, int nnrecv, int ncsend, int ncrecv, int *nshared, int *nsend, int *nrecv, int *csend, int *crecv)
Definition: FluRegion.H:286
vector< FluBorder > _borders
Definition: FluRegion.H:133
int _constr_type
Definition: FluRegion.H:108
bool ReadFluGridASCII(const std::string &p="", double=0.0, bool=true)
Definition: FluRegion.C:975
unsigned int _ntri
Definition: FluRegion.H:102
vector< int > _bflag
Definition: FluRegion.H:75
unsigned int _nnodes
Definition: FluRegion.H:126
void AddFluPatch(int patch_id, int ntri, int ngtri, int *tris, int nquad, int ngquad, int *quads)
Definition: FluRegion.H:307
double _resid
Definition: FluRegion.H:29
bool ReadFluCOM(const string &p="")
Definition: FluRegion.C:249
vector< double > _nf_alp
Definition: FluRegion.H:65
std::vector< unsigned int > _quadconn
Definition: GEM.H:132
string _casepath
Definition: FluRegion.H:122
unsigned int _nquad
Definition: FluRegion.H:103
unsigned int _cell_ordering[4]
Definition: GEM.H:303
unsigned int _npatches_total
Definition: FluRegion.H:131
unsigned int _ngpris
Definition: GEM.H:290
vector< double > _rhofvf_alp
Definition: FluRegion.H:64
bool RegisterSoln(const string &, bool)
Definition: FluRegion.C:1705
vector< double > _pf
Definition: FluRegion.H:67
vector< double > _Tflm_alp
Definition: FluRegion.H:73
vector< double > _qc
Definition: FluRegion.H:68
unsigned int nelem()
Definition: FluRegion.C:42
vector< double > _rhof
Definition: FluRegion.H:31
bool PopulateFluPatches(const string &pre, vector< DB > &db)
Definition: FluRegion.H:202
unsigned int _id
Definition: GEM.H:127
std::vector< unsigned int > _triconn
Definition: GEM.H:131
bool BuildPatchMapping(map< unsigned int, unsigned int > &, const string &p="")
Definition: FluRegion.C:1534
unsigned int _id
Definition: GEM.H:284
bool RegisterVolumeSoln(bool)
Definition: FluRegion.C:1689
std::ostream * _out
Definition: GEM.H:308
bool WriteGSPASCII(const std::string &p="", unsigned int=0, double=0.0, bool=true)
Definition: FluRegion.C:1389
std::ostream * _out
Definition: GEM.H:233
unsigned int _ngquad
Definition: GEM.H:129
void report()
Definition: FluRegion.C:47
unsigned int _nhex
Definition: FluRegion.H:128
bool _unsteady
Definition: FluRegion.H:30
bool WriteFluDIM(const string &p="", double=0.0, bool=true)
Definition: FluRegion.C:378
std::vector< unsigned int > _pyrconn
Definition: GEM.H:297
unsigned int _npris
Definition: FluRegion.H:130
bool validate_comm_list(int ncsend, int ncrecv, int *csend, int *crecv)
Definition: FluRegion.H:236
std::vector< unsigned int > _hexconn
Definition: GEM.H:301
void Resize(unsigned int ncells=0, unsigned int nfaces=0, unsigned int nnodes=0)
Definition: FluRegion.H:45
void populate(int rpid, int nnshared, int nnsend, int nnrecv, int ncsend, int ncrecv, int *sharedn, int *sendn, int *recvn, int *sendc, int *recvc)
Definition: GEM.C:193
bool ReadFluDIM(const string &p="", double=0.0, bool=true)
Definition: FluRegion.C:463
vector< pair< unsigned int, unsigned int > > _cellmap
Definition: FluRegion.H:135
bool WriteFluNative(const string &p="")
Definition: FluRegion.C:114
FluSurfSoln _soln
Definition: FluRegion.H:104
vector< FluPatch > _patches
Definition: FluRegion.H:132
vector< double > _mdot_alp
Definition: FluRegion.H:72
vector< double > _rhovf
Definition: FluRegion.H:32
vector< double > _tf
Definition: FluRegion.H:70
vector< double > _gsp
Definition: FluRegion.H:74
bool BuildCellMapping()
bool RegisterSurfaceSoln(bool)
Definition: FluRegion.C:1742
vector< double > _af
Definition: FluRegion.H:36
unsigned int _ntet
Definition: FluRegion.H:127
bool InitSurfaceSoln(const string &p="")
Definition: FluRegion.C:862
FluVolumeSoln _soln
Definition: FluRegion.H:134
FluSurfSoln(unsigned int nfaces=0, unsigned int nvert=0)
Definition: FluRegion.H:77
bool _debug
Definition: GEM.H:310
unsigned int _ngtet
Definition: GEM.H:287
bool Cell2Elem(unsigned int, unsigned int &, unsigned int &)
void Resize(unsigned int nfaces=0, unsigned int nvert=0)
Definition: FluRegion.H:81
double _current_time
Definition: FluRegion.H:28
bool InitVolumeSoln(const string &p="")
Definition: FluRegion.C:776
std::ostream * _out
Definition: GEM.H:140
vector< double > _disp
Definition: FluRegion.H:37
vector< double > _gsp
Definition: FluRegion.H:38
vector< double > _Tb_alp
Definition: FluRegion.H:71
FluRegion(const GEM_Partition &gp)
Definition: FluRegion.H:187
istream & SkipLines(istream &, unsigned int=1)
bool WriteFluGridASCII(const string &pre="", double=0.0, bool=true)
Definition: FluRegion.C:1085
std::vector< double > _nc
Definition: GEM.H:293
std::vector< unsigned int > _prisconn
Definition: GEM.H:299
vector< double > _rhof_alp
Definition: FluRegion.H:66
bool ReadGSPASCII(const std::string &p="", unsigned int=0, double=0.0, bool=true)
Definition: FluRegion.C:1480
unsigned int _ngpyr
Definition: GEM.H:289
unsigned int _nghex
Definition: GEM.H:288
bool RegisterFluSurfaceMesh()
Definition: FluRegion.C:1783
string _casename
Definition: FluRegion.H:121
FluVolumeSoln(unsigned int ncells=0, unsigned int nfaces=0, unsigned int nnodes=0)
Definition: FluRegion.H:39
unsigned int _ngtri
Definition: GEM.H:128
vector< double > _qr
Definition: FluRegion.H:69