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 490 491 492 493 494 495 496 497 498 499 500 501 502 503 504 505 506 507 508 509 510 511 512 513 514 515 516 517 518 519 520 521 522 523 524 525 526 527 528 529 530 531 532 533 534 535 536 537 538 539 540 541 542 543 544 545 546 547 548 549 550 551 552 553 554 555 556 557 558 559 560 561 562 563 564 565 566 567 568 569 570 571 572 573 574 575 576 577 578 579 580 581 582 583 584 585 586 587 588 589 590 591 592 593 594 595 596 597 598 599 600 601 602 603 604 605 606 607 608 609 610 611 612 613 614 615 616 617 618 619 620 621 622 623 624 625 626 627 628 629 630 631 632 633 634 635 636 637 638 639 640 641 642 643 644 645 646 647 648 649 650 651 652 653 654 655 656 657 658 659 660 661 662 663 664 665 666 667 668 669 670 671 672 673 674 675 676 677 678 679 680 681 682 683 684 685 686 687 688 689 690 691 692 693 694 695 696 697 698 699 700 701 702 703 704 705 706 707 708 709 710 711 712 713 714 715 716 717 718 719 720 721 722 723 724 725 726 727 728 729 730 731 732 733 734 735 736 737 738 739 740 741 742 743 744 745 746 747 748 749 750 751 752 753 754 755 756 757 758 759 760 761 762 763 764 765 766 767 768 769 770 771 772 773 774 775 776 777 778 779 780 781 782 783 784 785 786 787 788 789 790 791 792 793 794 795 796 797 798 799 800 801 802 803 804 805 806 807 808 809 810 811 812 813 814 815 816
|
#include <meshing.hpp>
#include "writeuser.hpp"
#ifdef NG_CGNS
#include <variant>
#include <cgnslib.h>
namespace netgen::cg
{
typedef ngcore::ClosedHashTable<ngcore::IVec<3,size_t>, size_t> PointTable;
int getDim(ElementType_t type)
{
switch(type)
{
case BAR_2:
case BAR_3:
return 1;
case TRI_3:
case TRI_6:
case QUAD_4:
case QUAD_8:
return 2;
case TETRA_4:
case TETRA_10:
case PYRA_5:
case PYRA_13:
case HEXA_8:
case HEXA_20:
case PENTA_6:
case PENTA_15:
return 3;
default:
throw Exception("Read CGNS: unknown element type " + string(cg_ElementTypeName(type)));
}
}
Segment ReadCGNSElement1D( ElementType_t type, FlatArray<cgsize_t> verts )
{
int np;
cg_npe(type, &np);
Segment s;
for (auto i : Range(np))
s[i] = verts[i];
return s;
}
Element2d ReadCGNSElement2D( ElementType_t type, FlatArray<cgsize_t> verts )
{
// static constexpr int map_tri3[] = {0,2,1};
static constexpr int map_tri6[] = {0,2,1,3,5,4}; // untested
// static constexpr int map_quad4[] = {0,3,2,1};
static constexpr int map_quad8[] = {0,3,2,1,4,7,6,5}; // untested
const int * map = nullptr;
switch(type)
{
case TRI_3:
// map = map_tri3;
break;
case QUAD_4:
// map = map_quad4;
break;
case TRI_6:
map = map_tri6;
break;
case QUAD_8:
map = map_quad8;
break;
default:
throw Exception("Read CGNS: unknown element type " + string(cg_ElementTypeName(type)));
}
int np;
cg_npe(type, &np);
Element2d el(np);
for (auto i : Range(np))
el[i] = verts[i];
return el;
}
Element ReadCGNSElement3D( ElementType_t type, FlatArray<cgsize_t> verts )
{
static constexpr int map_tet4[] = {0,2,1,3};
static constexpr int map_prism6[] = {0,2,1,3,5,4};
static constexpr int map_pyra5[] = {0,3,2,1,4};
static constexpr int map_hexa8[] = {0,3,2,1,4,7,6,5};
int np;
cg_npe(type, &np);
const int * map = nullptr;
switch(type)
{
case TETRA_4:
map = map_tet4; break;
case PYRA_5:
map = map_pyra5; break;
case PENTA_6:
map = map_prism6; break;
case HEXA_8:
map = map_hexa8; break;
// TODO: Second order elements
case TETRA_10:
case PYRA_13:
case HEXA_20:
case PENTA_15:
default:
throw Exception("Read CGNS: unknown element type " + string(cg_ElementTypeName(type)));
}
Element el(np);
for (auto i : Range(np))
el[i] = verts[map[i]];
return el;
}
void WriteCGNSElement( const Segment & el, Array<cgsize_t> & verts )
{
verts.Append(BAR_2);
verts.Append(el[0]);
verts.Append(el[1]);
}
void WriteCGNSElement( const Element2d & el, Array<cgsize_t> & verts )
{
static constexpr int map_tri6[] = {0,2,1,3,5,4}; // untested
static constexpr int map_quad8[] = {0,3,2,1,4,7,6,5}; // untested
ElementType_t type;
const int * map = nullptr;
switch(el.GetType())
{
case TRIG:
type = TRI_3;
break;
case QUAD:
type = QUAD_4;
break;
case TRIG6:
type = TRI_6;
map = map_tri6;
break;
case QUAD8:
type = QUAD_8;
map = map_quad8;
break;
// TODO: Second order elements
default:
throw Exception("Write CGNS: unknown element type " + ToString(el.GetType()));
}
verts.Append(type);
for (auto i : Range(el.GetNP()))
verts.Append(el[i]);
}
void WriteCGNSElement( const Element & el, Array<cgsize_t> & verts )
{
static constexpr int map_tet4[] = {0,2,1,3};
static constexpr int map_prism6[] = {0,2,1,3,5,4};
static constexpr int map_pyra5[] = {0,3,2,1,4};
static constexpr int map_hexa8[] = {0,3,2,1,4,7,6,5};
ElementType_t type;
const int * map = nullptr;
switch(el.GetType())
{
case TET:
map = map_tet4;
type = TETRA_4;
break;
case PYRAMID:
type = PYRA_5;
map = map_pyra5;
break;
case PRISM:
type = PENTA_6;
map = map_prism6;
break;
case HEX:
type = HEXA_8;
map = map_hexa8;
break;
// TODO: Second order elements
default:
throw Exception("Write CGNS: unknown element type " + ToString(el.GetType()));
}
verts.Append(type);
for (auto i : Range(el.GetNP()))
verts.Append(el[map[i]]);
}
int WriteCGNSRegion( const Mesh & mesh, int dim, int index, int fn, int base, int zone, int ne_before )
{
int meshdim = mesh.GetDimension();
int codim = meshdim-dim;
if(codim < 0 || codim > 2)
return 0;
// make sure that each material/boundary name is unique
string prefix[] = { "dom_", "bnd_", "bbnd_" };
string name = prefix[meshdim-dim] + ToString(index) + "_";
if(codim==0) name += mesh.GetMaterial(index+1);
if(codim==1) name += *mesh.GetBCNamePtr(index);
if(codim==2) name += mesh.GetCD2Name(index);
int ne = 0;
Array<int> data;
if(dim==3)
for(const auto el : mesh.VolumeElements())
if(el.GetIndex()==index)
{
ne++;
WriteCGNSElement(el, data);
}
if(dim==2)
for(const auto el : mesh.SurfaceElements())
if(el.GetIndex()==index)
{
ne++;
WriteCGNSElement(el, data);
}
if(dim==1)
for(const auto el : mesh.LineSegments())
if(el.si==index)
{
ne++;
WriteCGNSElement(el, data);
}
if(ne==0)
return 0;
int section;
// int start = 1;
// int end = ne;
#if CGNS_VERSION < 3400
cg_section_write(fn,base,zone, name.c_str(), MIXED, ne_before+1, ne_before+ne, 0, &data[0], §ion);
#else
cg_poly_section_write(fn,base,zone, name.c_str(), MIXED, ne_before+1, ne_before+ne, 0, &data[0], nullptr, §ion);
#endif
return ne;
}
// maps cgns node type to ngsolve node type
// enum NODE_TYPE { NT_VERTEX = 0, NT_EDGE = 1, NT_FACE = 2, NT_CELL = 3, NT_ELEMENT = 4, NT_FACET = 5 };
int getNodeType( GridLocation_t location )
{
switch(location)
{
case Vertex:
return 0;
case CellCenter:
return 3;
case FaceCenter:
return 2;
case EdgeCenter:
return 1;
default:
throw Exception("Read CGNS: unknown grid location " + string(cg_GridLocationName(location)));
}
}
GridLocation_t getCGNodeType( int node_type )
{
switch(node_type)
{
case 0:
return Vertex;
case 1:
return EdgeCenter;
case 2:
return FaceCenter;
case 3:
return CellCenter;
default:
throw Exception("Write CGNS: unknown node type " + ToString(node_type));
}
}
struct Solution
{
int fn, base, zone, solution;
string name;
GridLocation_t location; // solution is defined on either cells, faces, edges or vertices
PointSetType_t point_type;
cgsize_t n_points;
Array<string> field_names;
Array<DataType_t> field_datatypes;
Solution() = default;
Solution(int fn_, int base_, int zone_, int solution_)
: fn(fn_), base(base_), zone(zone_), solution(solution_)
{
char solname[100];
cg_sol_info(fn, base, zone, solution, solname, &location);
name = solname;
cg_sol_ptset_info(fn, base, zone, solution, &point_type, &n_points);
int n_fields = 0;
cg_nfields(fn, base, zone, solution, &n_fields);
field_names.SetSize(n_fields);
field_datatypes.SetSize(n_fields);
for(auto fi : Range(n_fields))
{
char buf[100];
cg_field_info(fn, base, zone, solution, fi+1, &field_datatypes[fi], buf);
field_names[fi] = buf;
}
}
};
struct Zone
{
ZoneType_t zone_type;
int fn, base, zone;
int first_index_1d, first_index_2d, first_index_3d;
int nv=0, ne_1d=0, ne_2d=0, ne_3d=0;
Array<string> names_1d, names_2d, names_3d;
string name;
cgsize_t size[3];
Array<Solution> solutions;
Zone(int fn_, int base_, int zone_)
: fn(fn_), base(base_), zone(zone_)
{
cg_zone_type(fn, base, zone, &zone_type);
char zone_name[100];
cg_zone_read(fn,base,zone, zone_name, size);
nv = size[0];
int n_solutions;
cg_nsols(fn, base, zone, &n_solutions);
solutions.SetSize(n_solutions);
for(auto si : Range(n_solutions))
solutions[si] = Solution{fn, base, zone, si+1};
}
void ReadSolutions( int meshdim, std::vector<string> & sol_names, std::vector<Array<double>> & sol_values, std::vector<int> & sol_locations )
{
static Timer tall("CGNS::ReadSolutions"); RegionTimer rtall(tall);
for (auto & sol : solutions)
{
for (auto fi : Range(sol.field_names.Size()))
{
cgsize_t size = sol.n_points;
size=0; // TODO: check if sol.point_type is a list or range, and handle appropriately
if(size==0)
{
switch(sol.location)
{
case Vertex:
size = nv;
break;
case CellCenter:
size = (meshdim == 3 ? ne_3d : ne_2d);
break;
case FaceCenter:
case IFaceCenter:
case JFaceCenter:
case KFaceCenter:
case EdgeCenter:
default:
throw Exception("Read CGNS: unknown grid location " + string(cg_GridLocationName(sol.location)));
}
}
auto values = Array<double>(size);
cgsize_t imin = 1UL;
cg_field_read(fn, base, zone, sol.solution, sol.field_names[fi].c_str(), RealDouble, &imin, &size, &values[0]);
sol_names.push_back(sol.field_names[fi]);
sol_values.emplace_back(std::move(values));
sol_locations.push_back(getNodeType(sol.location));
}
}
}
void ReadMesh( Mesh & mesh, PointTable & point_table )
{
static Timer tall("CGNS::ReadMesh-Zone"); RegionTimer rtall(tall);
static Timer tsection("CGNS::ReadMesh-Section");
first_index_1d = mesh.GetRegionNamesCD(2).Size();
first_index_2d = mesh.GetRegionNamesCD(1).Size();
first_index_3d = mesh.GetRegionNamesCD(0).Size();
Array<double> x(nv), y(nv), z(nv);
cgsize_t imin=1;
cg_coord_read(fn,base,zone, "CoordinateX", RealDouble, &imin, &nv, &x[0]);
cg_coord_read(fn,base,zone, "CoordinateY", RealDouble, &imin, &nv, &y[0]);
cg_coord_read(fn,base,zone, "CoordinateZ", RealDouble, &imin, &nv, &z[0]);
Array<cgsize_t> point_map(nv);
for(auto i : Range(nv))
{
ngcore::IVec<3,size_t> hash = {*reinterpret_cast<size_t*>(&x[i]), *reinterpret_cast<size_t*>(&y[i]), *reinterpret_cast<size_t*>(&z[i])};
size_t pi_ng;
size_t pos;
// check if this point is new
if( point_table.PositionCreate (hash, pos) )
{
pi_ng = mesh.AddPoint( {x[i], y[i], z[i]} );
point_table.SetData(pos, pi_ng);
}
else
point_table.GetData(pos, pi_ng);
point_map[i] = pi_ng;
}
int nsections;
cg_nsections(fn, base, zone, &nsections);
int index_1d = first_index_1d;
int index_2d = first_index_2d;
int index_3d = first_index_3d;
for (auto section : Range(1,nsections+1))
{
RegionTimer rtsection(tsection);
char sec_name[100];
ElementType_t type;
cgsize_t start, end;
int nbndry, parent_flag;
cg_section_read(fn, base, zone, section, sec_name, &type, &start, &end, &nbndry, &parent_flag);
PrintMessage(4, "Read section ", section, " with name ", sec_name, " and element type ", cg_ElementTypeName(type));
string ngname{sec_name};
for (char & c : ngname)
if(c==' ')
c = '_';
if(type==MIXED)
{
bool have_1d_elements = false;
bool have_2d_elements = false;
bool have_3d_elements = false;
cgsize_t nv;
cg_ElementDataSize(fn, base, zone, section, &nv);
Array<cgsize_t> vertices(nv);
#if CGNS_VERSION < 3400
cg_elements_read(fn, base, zone, section, &vertices[0], nullptr);
#else
cg_poly_elements_read(fn, base, zone, section, &vertices[0], nullptr, nullptr);
#endif
size_t vi = 0;
while(vi<nv)
{
auto type = static_cast<ElementType_t>(vertices[vi++]);
int dim = getDim(type);
int np;
cg_npe(type, &np);
for (auto & v : vertices.Range(vi, vi+np))
v = point_map[v-1];
if(dim==1)
{
if(!have_1d_elements)
{
index_1d++;
have_1d_elements = true;
mesh.AddEdgeDescriptor(EdgeDescriptor{});
names_1d.Append(ngname);
}
auto el = ReadCGNSElement1D(type, vertices.Range(vi, vertices.Size()));
el.si = index_1d;
mesh.AddSegment(el);
vi += el.GetNP();
ne_1d++;
}
if(dim==2)
{
if(!have_2d_elements)
{
index_2d++;
have_2d_elements = true;
mesh.AddFaceDescriptor(FaceDescriptor(index_2d, 1, 0, 1));
names_2d.Append(ngname);
}
auto el = ReadCGNSElement2D(type, vertices.Range(vi, vertices.Size()));
el.SetIndex(index_2d);
mesh.AddSurfaceElement(el);
vi += el.GetNP();
ne_2d++;
}
if(dim==3)
{
if(!have_3d_elements)
{
index_3d++;
have_3d_elements = true;
names_3d.Append(ngname);
}
auto el = ReadCGNSElement3D(type, vertices.Range(vi, vertices.Size()));
el.SetIndex(index_3d);
mesh.AddVolumeElement(el);
vi += el.GetNP();
ne_3d++;
}
}
}
else
{
int dim = getDim(type);
cgsize_t nv;
cg_ElementDataSize(fn, base, zone, section, &nv);
int np=0;
cg_npe(type, &np);
Array<cgsize_t> vertices(nv);
cg_elements_read(fn, base, zone, section, &vertices[0], nullptr);
for (auto & v : vertices)
v = point_map[v-1];
int ne_section = nv/np;
if(dim==1)
{
index_1d++;
mesh.AddEdgeDescriptor(EdgeDescriptor{});
names_1d.Append(ngname);
for(auto i : Range(ne_section))
{
auto el = ReadCGNSElement1D(type, vertices.Range(np*i, np*(i+1)));
el.si = index_1d;
mesh.AddSegment(el);
}
ne_1d += ne_section;
}
if(dim==2)
{
index_2d++;
mesh.AddFaceDescriptor(FaceDescriptor(index_2d, 1, 0, 1));
names_2d.Append(ngname);
for(auto i : Range(ne_section))
{
auto el = ReadCGNSElement2D(type, vertices.Range(np*i, np*(i+1)));
el.SetIndex(index_2d);
mesh.AddSurfaceElement(el);
}
ne_2d += ne_section;
}
if(dim==3)
{
index_3d++;
names_3d.Append(ngname);
for(auto i : Range(ne_section))
{
auto el = ReadCGNSElement3D(type, vertices.Range(np*i, np*(i+1)));
el.SetIndex(index_3d);
mesh.AddVolumeElement(el);
}
ne_3d += ne_section;
}
}
}
mesh.GetRegionNamesCD(2).SetSize(index_1d);
mesh.GetRegionNamesCD(1).SetSize(index_2d);
mesh.GetRegionNamesCD(0).SetSize(index_3d);
mesh.GetRegionNamesCD(2) = nullptr;
mesh.GetRegionNamesCD(1) = nullptr;
mesh.GetRegionNamesCD(0) = nullptr;
}
void SetNames( Mesh & mesh )
{
if(mesh.GetDimension() == 2)
{
for (auto i : Range(names_1d.Size()))
mesh.SetBCName(first_index_1d + i, names_1d[i]);
for (auto i : Range(names_2d.Size()))
mesh.SetMaterial(first_index_2d + i +1, names_2d[i]);
}
else
{
for (auto i : Range(names_1d.Size()))
mesh.SetCD2Name(first_index_1d + i +1, names_1d[i]);
for (auto i : Range(names_2d.Size()))
{
mesh.SetBCName(first_index_2d + i, names_2d[i]);
mesh.GetFaceDescriptor(first_index_2d + i +1).SetDomainIn(first_index_3d+1);
}
for (auto i : Range(names_3d.Size()))
mesh.SetMaterial(first_index_3d + i +1, names_3d[i]);
}
}
};
}
namespace netgen
{
int ReadCGNSMesh (Mesh & mesh, const filesystem::path & filename, Array<unique_ptr<cg::Zone>> & zones)
{
mesh.SetDimension(3);
static Timer tall("CGNS::ReadMesh"); RegionTimer rtall(tall);
int fn;
cg_open(filename.string().c_str(),CG_MODE_READ,&fn);
int base = 1;
int nzones;
cg_nzones(fn, base, &nzones);
int n_vertices = 0;
for (auto zi : Range(1, nzones+1))
{
int size[3];
char name[100];
cg_zone_read(fn,base,zi, name, size);
n_vertices += size[0];
}
cg::PointTable points(2*n_vertices);
for (auto zi : Range(1, nzones+1))
{
ZoneType_t zone_type;
cg_zone_type(fn, base, zi, &zone_type);
if(zone_type != Unstructured )
{
PrintMessage(2, "skipping zone with type ", cg_ZoneTypeName(zone_type) );
continue;
}
auto zone = make_unique<cg::Zone>(fn, base, zi);
zone->ReadMesh( mesh, points );
zones.Append(std::move(zone));
}
if(mesh.GetNE() == 0)
mesh.SetDimension(2);
for (auto & zone : zones)
zone->SetNames(mesh);
mesh.UpdateTopology();
const auto & topo = mesh.GetTopology();
for (auto sei : Range(mesh.SurfaceElements()))
{
int ei0, ei1;
topo.GetSurface2VolumeElement (sei+1, ei0, ei1);
auto si = mesh.SurfaceElement(sei).GetIndex();
auto & fd = mesh.GetFaceDescriptor(si);
if(ei0>0)
{
int i0 = mesh.VolumeElement(ei0).GetIndex();
if(fd.DomainIn()!=i0)
fd.SetDomainOut(i0);
}
if(ei1>0)
{
int i1 = mesh.VolumeElement(ei1).GetIndex();
if(fd.DomainIn()!=i1)
fd.SetDomainOut(i1);
}
}
return fn;
}
void ReadCGNSMesh (Mesh & mesh, const filesystem::path & filename)
{
Array<unique_ptr<cg::Zone>> zones;
int fn = ReadCGNSMesh(mesh, filename, zones);
cg_close(fn);
}
// Reads mesh and solutions of .csns file
tuple<shared_ptr<Mesh>, vector<string>, vector<Array<double>>, vector<int>> ReadCGNSFile(const filesystem::path & filename, int base)
{
static Timer tall("CGNS::ReadFile"); RegionTimer rtall(tall);
auto mesh = make_shared<Mesh>();
Array<unique_ptr<cg::Zone>> zones;
int fn = ReadCGNSMesh(*mesh, filename, zones);
std::vector<string> names;
std::vector<Array<double>> values;
std::vector<int> locations;
for (auto & zone : zones)
zone->ReadSolutions( mesh->GetDimension(), names, values, locations );
cg_close(fn);
return std::make_tuple(mesh, names, values, locations);
}
void WriteCGNSMesh (const Mesh & mesh, int fn, int & base, int & zone)
{
int dim = mesh.GetDimension();
cg_base_write(fn, "mesh", dim, dim, &base);
int nv = static_cast<int>(mesh.GetNV());
int ne = mesh.GetNE();
Array<double> x, y, z;
for(auto & p : mesh.Points())
{
x.Append(p[0]);
y.Append(p[1]);
z.Append(p[2]);
}
cgsize_t isize[3] = { nv, ne, 0 };
cg_zone_write(fn,base, "mesh", isize, Unstructured, &zone);
int coord;
cg_coord_write(fn,base,zone, RealDouble, "CoordinateX", &x[0], &coord);
cg_coord_write(fn,base,zone, RealDouble, "CoordinateY", &y[0], &coord);
cg_coord_write(fn,base,zone, RealDouble, "CoordinateZ", &z[0], &coord);
int imax3 = 0;
for(const auto & el : mesh.VolumeElements())
imax3 = max(imax3, el.GetIndex());
int imax2 = 0;
for(const auto & el : mesh.SurfaceElements())
imax2 = max(imax2, el.GetIndex());
int imax1 = 0;
for(const auto & el : mesh.LineSegments())
imax1 = max(imax1, el.si);
int ne_written = 0;
// int meshdim = mesh.GetDimension();
for(const auto i : IntRange(imax3))
ne_written += cg::WriteCGNSRegion(mesh, 3, i+1, fn, base, zone, ne_written);
for(const auto i : IntRange(imax2))
ne_written += cg::WriteCGNSRegion(mesh, 2, i+1, fn, base, zone, ne_written);
for(const auto i : IntRange(imax1))
ne_written += cg::WriteCGNSRegion(mesh, 1, i+1, fn, base, zone, ne_written);
}
void WriteCGNSMesh (const Mesh & mesh, const filesystem::path & filename)
{
static Timer tall("CGNS::WriteMesh"); RegionTimer rtall(tall);
int fn, base, zone;
cg_open(filename.string().c_str(),CG_MODE_WRITE,&fn);
WriteCGNSMesh(mesh, fn, base, zone);
cg_close(fn);
}
void WriteCGNSFile(shared_ptr<Mesh> mesh, const filesystem::path & filename, vector<string> fields, vector<Array<double>> values, vector<int> locations)
{
static Timer tall("CGNS::WriteFile"); RegionTimer rtall(tall);
int fn, base, zone;
cg_open(filename.string().c_str(),CG_MODE_WRITE,&fn);
WriteCGNSMesh(*mesh, fn, base, zone);
for(auto i : IntRange(fields.size()))
{
int section, field;
string name = "solution_" + ToString(i);
cg_sol_write(fn, base, zone, name.c_str(), cg::getCGNodeType(locations[i]), §ion);
cg_field_write(fn, base, zone, section, RealDouble, fields[i].c_str(), &values[i][0], &field);
}
cg_close(fn);
}
static RegisterUserFormat reg_cgns ("CGNS Format", {".cgns"},
static_cast<void(*)(Mesh &, const filesystem::path&)>(&ReadCGNSMesh),
static_cast<void(*)(const Mesh &, const filesystem::path&)>(&WriteCGNSMesh));
}
#endif // NG_CGNS
|