summaryrefslogtreecommitdiff
path: root/new/mesh.cpp
blob: 6180ea8023587f44f8eb31e7c86634226b3fd86f (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
// -*- mode: c++ -*-

#include "mesh.hpp"


// ------------------------------------------------------------------------
// -- PUBLIC METHODS
// ------------------------------------------------------------------------

void Mesh::build_surface_mesh(H5::H5File fp)
{
  std::cout << "[~>] Building the surface mesh" << std::endl;

  get_mesh_metadata(fp, SURFACE_NCOL);

  handle_surface_node_data(fp);

  build_surface_cells();
  
}

void Mesh::build_subsurface_mesh(H5::H5File fp)
{
  
  std::cout << "[~>] Building the subsurface mesh" << std::endl;

  get_mesh_metadata(fp, SUBSURFACE_NCOL);

  handle_subsurface_node_data(fp);

  build_subsurface_cells();
  
}

void Mesh::build_environment(H5::H5File ep, int flag)
{

  std::string groupname;
  std::string datasetname;
  
  // Read in velocities in x direction
    
  if (flag == SUBSURFACE_FLAG) {
    std::cout << "[~>] Building the subsurface environment" << std::endl;
    groupname = "darcy_velocity.cell.0";
    datasetname = SUBSURFACE_DATA_GROUP_NAME_H5;
  } else if (flag == SURFACE_FLAG) {
    std::cout << "[~>] Building the surface environment" << std::endl;
    groupname = "surface-velocity.cell.0";
    datasetname = SURFACE_DATA_GROUP_NAME_H5;
  }

  read_velocity(ep, groupname, datasetname, flag, 0);

  // Read in velocities in y direction

  if (flag == SUBSURFACE_FLAG) {
    groupname = "darcy_velocity.cell.1";
    datasetname = SUBSURFACE_DATA_GROUP_NAME_H5;
  } else if (flag == SURFACE_FLAG) {
    groupname = "surface-velocity.cell.1";
    datasetname = SURFACE_DATA_GROUP_NAME_H5;
  }

  read_velocity(ep, groupname, datasetname, flag, 1);

  // Read velocities in z direction
  
  if (flag == SUBSURFACE_FLAG) {
    groupname = "darcy_velocity.cell.2";
    datasetname = SUBSURFACE_DATA_GROUP_NAME_H5;
  } else if (flag == SURFACE_FLAG) {
    groupname = "surface-surface_subsurface_flux.cell.0";
    datasetname = SURFACE_DATA_GROUP_NAME_H5;
  }

  read_velocity(ep, groupname, datasetname, flag, 2);

  
  // Surface environment needs adjustments to the vertical velocity
  // and needs to read in water depth
  
  if (flag == SURFACE_FLAG) {
    
    // Rescale the exchange flux from mole to m/s
    for (int i = 0; i < n_cells_surface; i ++) {
      surface_cells[i].set_velocity(2, surface_cells[i].get_velocity(2) * 1.8e-5 * AREADIV);
      if (i % (n_cells_surface / MAXPRINT) == 0)
  	printf("[~~~~>] Surface velocity in cell[%d] = { %e, %e, %e } m/s\n", i,
	       surface_cells[i].get_velocity(0),
	       surface_cells[i].get_velocity(1),
	       surface_cells[i].get_velocity(2));
    }

    // Read in water depth (TODO)
    read_water_depth(ep);
    
  } else if (flag == SUBSURFACE_FLAG) {
    
    for (int i = 0; i < n_cells_subsurface; i ++) {
      if (i % (n_cells_subsurface / MAXPRINT) == 0)
  	printf("[~~~~>] Darcy velocity in cell[%d] = { %e, %e, %e } m/s\n", i,
	       subsurface_cells[i].get_velocity(0),
	       subsurface_cells[i].get_velocity(1),
	       subsurface_cells[i].get_velocity(2));
    }
    
  }
  
}

// ------------------------------------------------------------------------
// GETTERS AND SETTERS
// ------------------------------------------------------------------------

Mesh::Mesh()
{
  n_cells_surface = 0;
  n_cells_subsurface = 0;
  high_coordinates.resize(DIM);
  low_coordinates.resize(DIM);
}

int Mesh::get_n_cells(int flag)
{
  int n_cells;
  
  if (flag == SUBSURFACE_FLAG) {
    n_cells = n_cells_subsurface;
  } else {
    n_cells = n_cells_subsurface;
  }
  
  return n_cells;
}

Cell* Mesh::get_cell(int index, int flag)
{

  if (flag == SUBSURFACE_FLAG) {
    return &(subsurface_cells[index]);
  } else {
    return &(surface_cells[index]);
  } 
}

void Mesh::set_high_coordinate(int index, double value)
{
  high_coordinates[index] = value;
}

double Mesh::get_high_coordinate(int index)
{
  return high_coordinates[index];
}

void Mesh::set_low_coordinate(int index, double value)
{
  low_coordinates[index] = value;
}

double Mesh::get_low_coordinate(int index)
{
  return low_coordinates[index];
}



// ------------------------------------------------------------------------
// -- PRIVATE METHODS
// ------------------------------------------------------------------------

void Mesh::handle_surface_node_data(H5::H5File fp)
{

  H5::DataSet dataset_n = fp.openDataSet(MESH_NODES);
  H5::DataType datatype_n = dataset_n.getDataType();
  H5::DataSpace dataspace_n = dataset_n.getSpace();

  int rank = dataspace_n.getSimpleExtentNdims();

  hsize_t *dims_n = new hsize_t[rank];
  dataspace_n.getSimpleExtentDims(dims_n);

  for (int i = 0; i < rank; i ++)
    std::cout << "[~~~>] Dimension in " << i << " of nodes: " << dims_n[i] << std::endl;

  double *data_n = new double[dims_n[0]*dims_n[1]];
  dataset_n.read(data_n, H5::PredType::IEEE_F64LE);

  nodes.resize(dims_n[0]);

  for (size_t i = 0; i < dims_n[0]; i ++) {
    nodes[i].resize(dims_n[1]);
    for (size_t j = 0; j < dims_n[1]; j ++)
      nodes[i][j] = data_n[i * dims_n[1] + j];
  }
  
}

void Mesh::build_surface_cells()
{

  std::vector<double> coords1(DIM);
  std::vector<double> coords2(DIM);
  std::vector<double> coords3(DIM);
  std::vector<double> coords_avg(DIM);

  for (int i = 0; i < n_cells_surface; i ++) {

    int n1 = element_node_indices[i][1];
    int n2 = element_node_indices[i][2];
    int n3 = element_node_indices[i][3];

    std::vector<double> edge_length(DIM);

    for (size_t k = 0; k < DIM; k ++)
      edge_length[k] = DBL_MAX;

    for (size_t j = 0; j < DIM; j ++) {

      coords1[j] = nodes[n1][j];
      coords2[j] = nodes[n2][j];
      coords3[j] = nodes[n3][j];

      coords_avg[j] = (coords1[j] + coords2[j] + coords3[j]) / 3.0;

      edge_length[0] = fabs(coords1[j] - coords_avg[j]);
      edge_length[1] = fabs(coords2[j] - coords_avg[j]);
      edge_length[2] = fabs(coords3[j] - coords_avg[j]);
      
    }

    double smallest_distance = edge_length[0];
    double largest_distance  = edge_length[0];
    
    for (size_t k = 1; k < DIM; k ++) {
      smallest_distance = std::min(smallest_distance, edge_length[k]);
      largest_distance  = std::max(largest_distance, edge_length[k]);
    }

    surface_cells[i].set_characteristic_length(smallest_distance);
    surface_cells[i].set_maximum_length(largest_distance);

    for (size_t k = 0; k < DIM; k ++)
      surface_cells[i].set_coordinate(k, coords_avg[k]);

    surface_cells[i].set_cell_id((int) i);

    // find the coupled cell in the subsurface
    
    std::vector<double> distance(DIM);
    hsize_t coupled_neighbor = -1; // store current coupled neighbor here
    double dz = DBL_MAX;       // find the nearest cell in the vertical direction
    for (int k = 0; k < n_cells_subsurface; k ++) {

      for (size_t j = 0; j < DIM; j ++)
	distance[j] = subsurface_cells[k].get_coordinate(j) - coords_avg[j];

      double xy_distance = sqrt(distance[0]*distance[0] + distance[1]*distance[1]);

      if (xy_distance <= FUZZ) {
	dz = std::min(dz, fabs(distance[2]));
	coupled_neighbor = (hsize_t) k;
	if (dz <= 2.0 * subsurface_cells[k].get_characteristic_length())
	  break;
      }
      
    }

    surface_cells[i].set_coupled_neighbor(coupled_neighbor);
    subsurface_cells[coupled_neighbor].set_coupled_neighbor(i);

    if (i%(n_cells_surface/5)==0) {
      std::cout << "[~~~>] Coupling surface and subsurface meshes" << std::endl;
      std::cout << "[~~~~>] Subsurface neighbor of cell " << i << " is cell " << surface_cells[i].get_coupled_neighbor() << std::endl;
      std::cout << "[~~~~>] Surface neighbor of cell " << coupled_neighbor << " is cell " << subsurface_cells[coupled_neighbor].get_coupled_neighbor() << std::endl;
    }
    
  }

}

void Mesh::handle_subsurface_node_data(H5::H5File fp)
{
  
  H5::DataSet dataset_n = fp.openDataSet(MESH_NODES);
  H5::DataType datatype_n = dataset_n.getDataType();
  H5::DataSpace dataspace_n = dataset_n.getSpace();

  int rank = dataspace_n.getSimpleExtentNdims();

  hsize_t *dims_n = new hsize_t[rank];
  dataspace_n.getSimpleExtentDims(dims_n);

  for (int i = 0; i < rank; i ++)
    std::cout << "[~~~>] Dimension in " << i << " of nodes: " << dims_n[i] << std::endl;

  double *data_n = new double[dims_n[0]*dims_n[1]];
  dataset_n.read(data_n, H5::PredType::IEEE_F64LE);

  nodes.resize(dims_n[0]);

  for (size_t i = 0; i < dims_n[0]; i ++) {
    nodes[i].resize(dims_n[1]);
    for (size_t j = 0; j < dims_n[1]; j ++)
      nodes[i][j] = data_n[i * dims_n[1] + j];
  }
  
}

void Mesh::build_subsurface_cells()
{

  std::vector<double> coords1(DIM);
  std::vector<double> coords2(DIM);
  std::vector<double> coords3(DIM);
  std::vector<double> coords4(DIM);
  std::vector<double> coords5(DIM);
  std::vector<double> coords6(DIM);

  std::vector<double> coords_avg(DIM);

  for (size_t i = 0; i < DIM; i ++) {
    set_high_coordinate(i, DBL_MIN);
    set_low_coordinate(i, DBL_MAX);
  }

  for (int i = 0; i < n_cells_subsurface; i ++) {

    int n1 = element_node_indices[i][1];
    int n2 = element_node_indices[i][2];
    int n3 = element_node_indices[i][3];
    int n4 = element_node_indices[i][4];
    int n5 = element_node_indices[i][5];
    int n6 = element_node_indices[i][6];

    std::vector<double> edge_length(N_VERTICES);

    for (size_t k = 0; k < N_VERTICES; k ++)
      edge_length[k] = DBL_MAX;
    
    for (size_t j = 0; j < DIM; j++) {
      
      coords1[j] = nodes[n1][j];
      coords2[j] = nodes[n2][j];
      coords3[j] = nodes[n3][j];
      coords4[j] = nodes[n4][j];
      coords5[j] = nodes[n5][j];
      coords6[j] = nodes[n6][j];

      coords_avg[j] = (coords1[j] + coords2[j] + coords3[j] + coords4[j] + coords5[j] + coords6[j]) / 6.0;

      // the characteristic length is the radius of the insphere of
      // the wedge
      edge_length[0] = fabs(coords1[j] - coords_avg[j]);
      edge_length[1] = fabs(coords2[j] - coords_avg[j]);
      edge_length[2] = fabs(coords3[j] - coords_avg[j]);
      edge_length[3] = fabs(coords4[j] - coords_avg[j]);
      edge_length[4] = fabs(coords5[j] - coords_avg[j]);
      edge_length[5] = fabs(coords6[j] - coords_avg[j]);
      
    }

    double smallest_distance = edge_length[0];
    double largest_distance  = edge_length[0];
    
    for (size_t k = 1; k < N_VERTICES; k ++) {
      smallest_distance = std::min(smallest_distance, edge_length[k]);
      largest_distance  = std::max(largest_distance, edge_length[k]);
    }

    subsurface_cells[i].set_characteristic_length(smallest_distance);
    subsurface_cells[i].set_maximum_length(largest_distance);
    
    for (size_t k = 0; k < DIM; k ++) {

      // update lowest and highest corner of mesh
      high_coordinates[k] = std::max(coords_avg[k], high_coordinates[k]);
      low_coordinates[k] = std::min(coords_avg[k], low_coordinates[k]);

      subsurface_cells[i].set_coordinate(k, coords_avg[k]);
      
    }
    
    subsurface_cells[i].set_cell_id((int) i);
    
  }
  
}

void Mesh::get_mesh_metadata(H5::H5File fp, int flag)
{

  H5::DataSet dataset_e = fp.openDataSet(MESH_MIXED_ELEMENTS);
  H5::DataType datatype_e = dataset_e.getDataType();
  H5::DataSpace dataspace_e = dataset_e.getSpace();

  int rank = dataspace_e.getSimpleExtentNdims();

  hsize_t *dims_e = new hsize_t[rank];
  dataspace_e.getSimpleExtentDims(dims_e);

  for (int i = 0; i < rank; i ++)
      std::cout << "[~~~>] Dimension in " << i << " of element node map: " << dims_e[i] << std::endl;

  int *data_e = new int[dims_e[0]*dims_e[1]];
  dataset_e.read(data_e, H5::PredType::STD_I32LE);

  size_t nx_e = flag;
  size_t ny_e = dims_e[0] * dims_e[1] / nx_e;
  
  element_node_indices.resize(ny_e);

  for (size_t i = 0; i < ny_e; i ++) {
    element_node_indices[i].resize(nx_e);
    for (size_t j = 0; j < nx_e; j ++)
      element_node_indices[i][j] = data_e[i * nx_e + j];
  }
  
  if (flag == SUBSURFACE_NCOL) {
    n_cells_subsurface = ny_e;
    subsurface_cells.resize(ny_e);
  } else if (flag == SURFACE_NCOL) {
    n_cells_surface = ny_e;
    surface_cells.resize(ny_e);
  }

  std::cout << "[~~~>] Number of cells: " << ny_e << std::endl;
  
}

void Mesh::read_velocity(H5::H5File fp, std::string groupname, std::string datasetname, int flag, int index)
{
  
  H5::Group group		= fp.openGroup(groupname.c_str());
  H5::DataSet dataset		= group.openDataSet(datasetname.c_str());
  H5::DataType datatype		= dataset.getDataType();
  H5::DataSpace dataspace	= dataset.getSpace();
  
  int rank = dataspace.getSimpleExtentNdims();

  hsize_t *dims = new hsize_t[rank];
  dataspace.getSimpleExtentDims(dims);

  for (int i = 0; i < rank; i ++)
    std::cout << "[~~~>] Dimensions in " << i << " of velocity in x: " << dims[i] << std::endl;
  
  double *data = new double[dims[0]*dims[1]];
  dataset.read(data, H5::PredType::IEEE_F64LE);

  if (flag == SUBSURFACE_FLAG) {
    for (int i = 0; i < n_cells_subsurface; i ++)
      subsurface_cells[i].set_velocity(index, data[i]);
  } else if (flag == SURFACE_FLAG) {
    for (int i = 0; i < n_cells_surface; i ++)
      surface_cells[i].set_velocity(index, data[i]);
  }

}

void Mesh::read_water_depth(H5::H5File fp)
{

  H5::Group group		= fp.openGroup("surface-ponded_depth.cell.0");
  H5::DataSet dataset		= group.openDataSet(SURFACE_DATA_GROUP_NAME_H5);
  H5::DataType datatype		= dataset.getDataType();
  H5::DataSpace dataspace	= dataset.getSpace();
  
  int rank = dataspace.getSimpleExtentNdims();

  hsize_t *dims = new hsize_t[rank];
  dataspace.getSimpleExtentDims(dims);

  for (int i = 0; i < rank; i ++)
    std::cout << "[~~~>] Dimensions in " << i << " of velocity in x: " << dims[i] << std::endl;
  
  double *data = new double[dims[0]*dims[1]];
  dataset.read(data, H5::PredType::IEEE_F64LE);

  for (int i = 0; i < n_cells_surface; i ++) {
    surface_cells[i].set_water_depth(data[i]);
    if (i%(n_cells_surface/5)==0) {
      printf("[~~~~>] Surface ponded depth in cell[%d] =  %e m\n", i,
	     surface_cells[i].get_water_depth());
    }
  }
  
}