39 #define PETSC_HAVE_BROKEN_RECURSIVE_MACRO 41 #include "AbstractTetrahedralMeshWriter.hpp" 42 #include "AbstractTetrahedralMesh.hpp" 43 #include "DistributedTetrahedralMesh.hpp" 44 #include "MixedDimensionMesh.hpp" 45 #include "Version.hpp" 50 const char*
MeshEventHandler::EventName[] = {
"Tri write",
"BinTri write",
"VTK write",
"PVTK write",
"node write",
"ele write",
"face write",
"ncl write",
"comm1",
"comm2",
"Total"};
55 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
70 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
72 const std::string& rBaseName,
73 const bool clearOutputDir)
76 mNodesPerElement(ELEMENT_DIM+1),
77 mNodesPerBoundaryElement(ELEMENT_DIM),
79 mpDistributedMesh(nullptr),
82 mNodeCounterForParallelMesh(0),
83 mElementCounterForParallelMesh(0),
84 mBoundaryElementCounterForParallelMesh(0),
85 mCableElementCounterForParallelMesh(0),
86 mFilesAreBinary(false)
90 mpIters->pBoundaryElemIter =
nullptr;
93 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
104 if (
mpIters->pBoundaryElemIter)
106 delete mpIters->pBoundaryElemIter;
117 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
124 std::vector<double> coords(SPACE_DIM);
131 for (
unsigned j=0; j<SPACE_DIM; j++)
133 coords[j] = (*(
mpIters->pNodeIter))->GetPoint()[j];
148 status.MPI_ERROR = MPI_SUCCESS;
150 boost::scoped_array<double> raw_coords(
new double[SPACE_DIM]);
152 assert(status.MPI_ERROR == MPI_SUCCESS);
153 for (
unsigned j=0; j<coords.size(); j++)
155 coords[j] = raw_coords[j];
158 mNodeCounterForParallelMesh++;
167 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
182 for (
unsigned j=0; j<elem_data.
NodeIndices.size(); j++)
184 unsigned old_index = (*(
mpIters->pElemIter))->GetNodeGlobalIndex(j);
227 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
242 for (
unsigned j=0; j<boundary_elem_data.
NodeIndices.size(); j++)
244 unsigned old_index = (*(*(
mpIters->pBoundaryElemIter)))->GetNodeGlobalIndex(j);
249 ++(*(
mpIters->pBoundaryElemIter));
250 return boundary_elem_data;
259 assert(boundary_elem_data.
NodeIndices.size() == ELEMENT_DIM);
260 assert(!p_boundary_element->IsDeleted());
263 for (
unsigned j=0; j<ELEMENT_DIM; j++)
265 boundary_elem_data.
NodeIndices[j] = p_boundary_element->GetNodeGlobalIndex(j);
267 boundary_elem_data.
AttributeValue = p_boundary_element->GetAttribute();
277 return boundary_elem_data;
286 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
308 for (
unsigned j=0; j<2; j++)
330 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
333 bool invertMeshPermutation)
336 unsigned max_elements_all;
344 MPI_Allreduce(&max_elements_per_process, &max_elements_all, 1, MPI_UNSIGNED, MPI_MAX, PETSC_COMM_WORLD);
347 std::string node_connect_list_file_name = this->
mBaseName +
".ncl";
350 node_connect_list_file_name +=
"-tmp";
355 out_stream p_ncl_file = out_stream(
nullptr);
364 *p_ncl_file << max_elements_all <<
"\t";
365 *p_ncl_file <<
"\tBIN\n";
374 unsigned default_marker = UINT_MAX;
382 std::set<unsigned>& r_elem_set = iter->rGetContainingElementIndices();
383 std::vector<unsigned> elem_vector(r_elem_set.begin(), r_elem_set.end());
384 std::sort(elem_vector.begin(), elem_vector.end());
386 for (
unsigned elem_index=elem_vector.size(); elem_index<max_elements_all; elem_index++)
388 elem_vector.push_back(default_marker);
390 assert(elem_vector.size() == max_elements_all);
392 if (max_elements_all > 0u)
394 p_ncl_file->write((
char*)&elem_vector[0], elem_vector.size()*
sizeof(
unsigned));
410 const std::string real_node_connect_list_file_name = this->
mBaseName +
".ncl";
413 std::ifstream temp_ncl_file(temp_ncl_path.
GetAbsolutePath().c_str(), std::ios::binary);
415 std::string header_line;
416 getline(temp_ncl_file, header_line,
'\n');
417 (*p_ncl_file) << header_line <<
"\n";
418 const std::streampos data_start = temp_ncl_file.tellg();
419 const std::streamoff item_width = max_elements_all *
sizeof(
unsigned);
421 std::vector<unsigned> elem_vector(max_elements_all);
422 for (
unsigned node_index=0; node_index<rMesh.
GetNumAllNodes(); node_index++)
425 temp_ncl_file.seekg(data_start + item_width * permuted_index, std::ios_base::beg);
426 temp_ncl_file.read((
char*)&elem_vector[0], max_elements_all*
sizeof(
unsigned));
427 p_ncl_file->write((
char*)&elem_vector[0], max_elements_all*
sizeof(
unsigned));
440 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
443 bool keepOriginalElementIndexing)
454 mpIters->pNodeIter =
new NodeIterType(
mpMesh->GetNodeIteratorBegin());
457 mpIters->pElemIter =
new ElemIterType(
mpMesh->GetElementIteratorBegin());
460 mpIters->pBoundaryElemIter =
new BoundaryElemIterType(
mpMesh->GetBoundaryElementIteratorBegin());
463 if ((*(
mpIters->pElemIter)) !=
mpMesh->GetElementIteratorEnd())
469 if ((*(
mpIters->pBoundaryElemIter)) !=
mpMesh->GetBoundaryElementIteratorEnd())
502 unsigned node_map_index = 0;
503 if (
mpMesh->IsMeshChanging())
506 for (NodeIterType it =
mpMesh->GetNodeIteratorBegin(); it !=
mpMesh->GetNodeIteratorEnd(); ++it)
518 delete mpIters->pBoundaryElemIter;
519 mpIters->pBoundaryElemIter =
nullptr;
522 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
531 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
534 if (keepOriginalElementIndexing)
545 boost::scoped_array<double> raw_coords(
new double[SPACE_DIM]);
548 for (NodeIterType it =
mpMesh->GetNodeIteratorBegin(); it !=
mpMesh->GetNodeIteratorEnd(); ++it)
550 for (
unsigned j=0; j<SPACE_DIM; j++)
552 raw_coords[j] = it->GetPoint()[j];
554 MPI_Ssend(raw_coords.get(), SPACE_DIM, MPI_DOUBLE, 0, it->GetIndex(), PETSC_COMM_WORLD);
563 for (ElementIterType it =
mpMesh->GetElementIteratorBegin(); it !=
mpMesh->GetElementIteratorEnd(); ++it)
565 unsigned index = it->GetIndex();
570 raw_indices[j] = it->GetNodeGlobalIndex(j);
579 if (ELEMENT_DIM != 1)
581 boost::scoped_array<unsigned> raw_face_indices(
new unsigned[ELEMENT_DIM]);
583 for (BoundaryElementIterType it =
mpMesh->GetBoundaryElementIteratorBegin(); it !=
mpMesh->GetBoundaryElementIteratorEnd(); ++it)
585 unsigned index = (*it)->GetIndex();
586 if (
mpDistributedMesh->CalculateDesignatedOwnershipOfBoundaryElement(index) ==
true)
588 for (
unsigned j=0; j<ELEMENT_DIM; j++)
590 raw_face_indices[j] = (*it)->GetNodeGlobalIndex(j);
603 for (CableElementIterType it =
mpMixedMesh->GetCableElementIteratorBegin(); it !=
mpMixedMesh->GetCableElementIteratorEnd(); ++it)
605 unsigned index =(*it)->GetIndex();
606 if (
mpMixedMesh->CalculateDesignatedOwnershipOfCableElement(index) ==
true)
608 unsigned raw_cable_element_indices[2];
609 for (
unsigned j=0; j<2; j++)
611 raw_cable_element_indices[j] = (*it)->GetNodeGlobalIndex(j);
645 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
656 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
667 template<
unsigned ELEMENT_DIM,
unsigned SPACE_DIM>
virtual ElementData GetNextCableElement()
ElementData GetNextCableElement()
unsigned GetNodeGlobalIndex(unsigned localIndex) const
void UnpackElement(ElementData &rElementData, unsigned globalIndex, unsigned numIndices, unsigned tag)
const std::vector< unsigned > & rGetNodePermutation() const
std::string GetAbsolutePath() const
ElementData GetNextBoundaryElement()
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpMesh
static void BeginEvent(unsigned event)
AbstractTetrahedralMeshWriter(const std::string &rDirectory, const std::string &rBaseName, const bool clearOutputDir=true)
unsigned mNumBoundaryElements
void WriteFilesUsingMeshReader(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader)
NodeIterator GetNodeIteratorEnd()
virtual void CreateFilesWithHeaders()
unsigned CalculateMaximumContainingElementsPerProcess() const
virtual unsigned GetNumNodes()
OutputFileHandler * mpOutputFileHandler
virtual unsigned GetNumNodes() const
AbstractMesh< ELEMENT_DIM, SPACE_DIM >::NodeIterator * pNodeIter
unsigned mNodesPerBoundaryElement
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::BoundaryElementIterator * pBoundaryElemIter
virtual unsigned GetNumAllNodes() const
virtual std::vector< double > GetNextNode()
AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM >::ElementIterator * pElemIter
MeshWriterIterators< ELEMENT_DIM, SPACE_DIM > * mpIters
std::vector< unsigned > NodeIndices
void SetNewIndex(unsigned oldIndex, unsigned newIndex)
out_stream OpenOutputFile(const std::string &rFileName, std::ios_base::openmode mode=std::ios::out|std::ios::trunc) const
void WriteNclFile(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool invertMeshPermutation=false)
std::vector< double > GetNextNode()
static const char * EventName[11]
FileFinder FindFile(std::string leafName) const
DistributedTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > * mpDistributedMesh
virtual void WriteFilesUsingParallelMesh(bool keepOriginalElementIndexing=true)
ElementData GetNextElement()
virtual ~AbstractTetrahedralMeshWriter()
void PostElement(unsigned globalIndex, unsigned indices[], unsigned numIndices, unsigned tag, double attribute)
NodeIterator GetNodeIteratorBegin(bool skipDeletedNodes=true)
unsigned mElementCounterForParallelMesh
virtual void WriteFilesFooter()
MixedDimensionMesh< ELEMENT_DIM, SPACE_DIM > * mpMixedMesh
virtual void WriteFiles()=0
unsigned mNodesPerElement
unsigned GetNewIndex(unsigned oldIndex) const
virtual ElementData GetNextBoundaryElement()
unsigned mNodeCounterForParallelMesh
static std::string GetProvenanceString()
static void EndEvent(unsigned event)
std::vector< BoundaryElement< ELEMENT_DIM-1, SPACE_DIM > * >::const_iterator BoundaryElementIterator
unsigned mNumCableElements
virtual ElementData GetNextElement()
virtual void WriteFilesUsingMesh(AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh, bool keepOriginalElementIndexing=true)
unsigned mBoundaryElementCounterForParallelMesh
unsigned mCableElementCounterForParallelMesh
virtual void AppendLocalDataToFiles()
void WriteFilesUsingMeshReaderAndMesh(AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > &rMeshReader, AbstractTetrahedralMesh< ELEMENT_DIM, SPACE_DIM > &rMesh)
std::vector< Element< 1, SPACE_DIM > * >::const_iterator CableElementIterator
AbstractMeshReader< ELEMENT_DIM, SPACE_DIM > * mpMeshReader