Loading...
Searching...
No Matches
POD_example_00002.cpp

POD basis computation using voloctree. This example computes the POD basis starting from a database of simulations defined on the same mesh and evaluate the reconstruction of a snapshot not included in the database testing the hybrid interface. The reconstruction is performed by methods based on POD field structures and on PiercedStorage containers. Finally, the reconstruction of a PiercedStorage field is performed by using an adapted mesh field different from the POD mesh. To run: ./POD_example_00002
.

POD basis computation using voloctree. This example computes the POD basis starting from a database of simulations defined on the same mesh and evaluate the reconstruction of a snapshot not included in the database testing the hybrid interface. The reconstruction is performed by methods based on POD field structures and on PiercedStorage containers. Finally, the reconstruction of a PiercedStorage field is performed by using an adapted mesh field different from the POD mesh. To run: ./POD_example_00002

/*---------------------------------------------------------------------------*\
*
* bitpit
*
* Copyright (C) 2015-2021 OPTIMAD engineering Srl
*
* -------------------------------------------------------------------------
* License
* This file is part of bitpit.
*
* bitpit is free software: you can redistribute it and/or modify it
* under the terms of the GNU Lesser General Public License v3 (LGPL)
* as published by the Free Software Foundation.
*
* bitpit is distributed in the hope that it will be useful, but WITHOUT
* ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
* FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
* License for more details.
*
* You should have received a copy of the GNU Lesser General Public License
* along with bitpit. If not, see <http://www.gnu.org/licenses/>.
*
\*---------------------------------------------------------------------------*/
#include <array>
#if BITPIT_ENABLE_MPI
#include <mpi.h>
#endif
#include "pod.hpp"
#include "bitpit_patchkernel.hpp"
#include "bitpit_voloctree.hpp"
using namespace bitpit;
void run(int rank, int nProcs)
{
if (rank==0)
std::cout<< "0. Compute POD ..." << std::endl;
POD pod;
for (int i=1; i<11; i++)
pod.addSnapshot("./data", "test."+std::to_string(i));
pod.setMeshType(POD::MeshType::VOLOCTREE);
pod.setStaticMesh(true);
pod.setWriteMode(POD::WriteMode::DEBUG);
pod.setMemoryMode(POD::MemoryMode::MEMORY_NORMAL);
pod.setModeCount(9);
pod.setDirectory("pod");
pod.setName("pod.test.solver");
pod.run();
if (rank==0)
std::cout<< "1. Restore POD ..." << std::endl;
POD podr;
podr.setDirectory("pod");
podr.setName("pod.test.solver");
podr.restore();
if (rank==0)
std::cout<< "2. Read snapshot for reconstruction ..." << std::endl;
#if BITPIT_ENABLE_MPI
VolumeKernel* meshr = new VolOctree(MPI_COMM_NULL);
#else
VolumeKernel* meshr = new VolOctree();
#endif
{
int dumpBlock = (nProcs > 1) ? rank : -1;
std::string filename = "./data/test.0.mesh";
IBinaryArchive binaryReader(filename, dumpBlock);
meshr->restore(binaryReader.getStream());
binaryReader.close();
}
pod::PODField fieldr, reconr;
std::size_t nf;
int dumpBlock = (nProcs > 1) ? rank : -1;
std::string filename = "./data/test.0.data";
IBinaryArchive binaryReader(filename, dumpBlock);
std::istream &dataStream = binaryReader.getStream();
fieldr.mask = std::unique_ptr<PiercedStorage<bool>>(new PiercedStorage<bool>(1, &meshr->getCells()));
fieldr.mask->restore(dataStream);
std::size_t nsf;
utils::binary::read(dataStream, nsf);
std::vector<std::string> namesf;
namesf.resize(nsf);
for (std::size_t i = 0; i <nsf; ++i){
utils::binary::read(dataStream, namesf[i]);
}
fieldr.scalar = std::unique_ptr<pod::ScalarStorage>(new pod::ScalarStorage(nsf, &meshr->getCells()));
fieldr.scalar->restore(dataStream);
std::size_t nvf;
utils::binary::read(dataStream, nvf);
std::vector<std::array<std::string,3>> namevf;
namevf.resize(nvf);
for (std::size_t i = 0; i < nvf; ++i){
utils::binary::read(dataStream, namevf[i]);
}
fieldr.vector = std::unique_ptr<pod::VectorStorage>(new pod::VectorStorage(nvf, &meshr->getCells()));
fieldr.vector->restore(dataStream);
binaryReader.close();
nf=nsf+nvf;
pod.setSensorMask(*(fieldr.mask));
pod.reconstructFields(fieldr,reconr);
if (rank==0){
std::cout<< ">> N modes:" << pod.getModeCount() << std::endl;
std::cout<< ">> Reconstruction coeffs:" << std::endl;
for (std::size_t i = 0; i < nf; ++i)
std::cout<< pod.getReconstructionCoeffs()[i] << std::endl;
}
PiercedStorage<double> gfield0(nsf+3*nvf+1, &meshr->getCells(), PiercedSyncMaster::SyncMode::SYNC_MODE_JOURNALED);
PiercedStorage<double> gfieldAdapt(nsf+3*nvf+1, &meshr->getCells(), PiercedSyncMaster::SyncMode::SYNC_MODE_JOURNALED);
std::unordered_set<long> targetCells;
std::map<std::string, std::size_t> fields;
for (Cell & cell : meshr->getCells()){
long id = cell.getId();
for (std::size_t ifield=0; ifield<nsf; ifield++){
gfield0.at(id, ifield) = fieldr.scalar->at(id, ifield);
gfieldAdapt.at(id, ifield) = fieldr.scalar->at(id, ifield);
}
for (std::size_t ifield=0; ifield<nvf; ifield++){
for (std::size_t j=0; j<3; j++){
gfield0.at(id, ifield*3+nsf+j) = fieldr.vector->at(id, ifield)[j];
gfieldAdapt.at(id, ifield*3+nsf+j) = fieldr.vector->at(id, ifield)[j];
}
}
if (fieldr.mask->at(id)){
targetCells.insert(id);
gfield0.at(id, nsf+3*nvf) = fieldr.mask->at(id);
gfieldAdapt.at(id, nsf+3*nvf) = fieldr.mask->at(id);
}
}
for (std::size_t ifield=0; ifield<nsf; ifield++){
fields[namesf[ifield]] = ifield;
}
for (std::size_t ifield=0; ifield<nvf; ifield++){
for (std::size_t j=0; j<3; j++){
fields[namevf[ifield][j]] = ifield*3+nsf+j;
}
}
pod.reconstructFields(gfield0, meshr, fields, &targetCells);
if (rank==0){
std::cout<< " " << std::endl;
std::cout<< ">> Reconstruction coeffs hybrid interface:" << std::endl;
for (std::size_t i = 0; i < nf; ++i)
std::cout<< pod.getReconstructionCoeffs()[i] << std::endl;
std::cout<< " " << std::endl;
std::cout<< ">> Reconstruction error hybrid interface" << std::endl;
double maxerr = 0.0;
for (Cell & cell : meshr->getCells()){
long id = cell.getId();
for (std::size_t i = 0; i < nsf; ++i)
maxerr = std::max(maxerr, std::abs(gfield0.at(id, i) - fieldr.scalar->at(id, i)));
for (std::size_t i = 0; i < nvf; ++i){
for (std::size_t j=0; j<3; j++)
maxerr = std::max(maxerr, std::abs(gfield0.at(id, nsf+3*i+j) - fieldr.vector->at(id, i)[j]));
}
}
std::cout<< ">> Max error :" << maxerr << std::endl;
}
long nCells = meshr->getCellCount();
log::cout() << std::endl;
log::cout() << ">> Marking the cells to adapt... " << std::endl;
for (int i = 0; i < 100; ++i) {
long cellId = rand() % nCells * 2;
if (!meshr->getCells().exists(cellId)) {
continue;
}
for (auto neighId : meshr->findCellNeighs(cellId)) {
meshr->markCellForRefinement(neighId);
}
}
for (int i = 0; i < 50; ++i) {
long cellId = rand() % nCells * 2;
if (!meshr->getCells().exists(cellId)) {
continue;
}
if (fieldr.mask->at(cellId)){
meshr->markCellForCoarsening(cellId);
for (auto neighId : meshr->findCellNeighs(cellId)) {
if (fieldr.mask->at(neighId))
meshr->markCellForCoarsening(neighId);
}
}
}
log::cout() << std::endl;
log::cout() << ">> Initial number of cells... " << nCells << std::endl;
{
std::vector<adaption::Info> preadaptInfo = meshr->adaptionPrepare(true);
for (adaption::Info & info : preadaptInfo){
for (long & id : info.previous){
oldData.insert(id, std::vector<double>(nsf+3*nvf+1, 0.0));
for (std::size_t i=0; i<nsf+3*nvf+1; i++)
oldData[id][i] = gfield0.at(id,i);
}
}
std::vector<adaption::Info> adaptInfo = meshr->adaptionAlter(true);
nCells = meshr->getCellCount();
log::cout() << ">> Final number of cells... " << nCells << std::endl;
for (adaption::Info & info : adaptInfo){
if (info.type == adaption::Type::TYPE_RENUMBERING){
for (std::size_t i=0; i<nsf+3*nvf+1; i++){
gfield0.at(info.current[0],i) = oldData[info.previous[0]][i];
gfieldAdapt.at(info.current[0],i) = oldData[info.previous[0]][i];
}
}
if (info.type == adaption::Type::TYPE_REFINEMENT){
for (long & id : info.current){
for (std::size_t i=0; i<nsf+3*nvf+1; i++){
gfield0.at(id,i) = oldData[info.previous[0]][i];
gfieldAdapt.at(id,i) = oldData[info.previous[0]][i];
}
}
}
if (info.type == adaption::Type::TYPE_COARSENING){
for (std::size_t i=0; i<nsf+3*nvf+1; i++){
gfield0.at(info.current[0], i) = 0.0;
gfieldAdapt.at(info.current[0], i) = 0.0;
}
for (long & id : info.previous){
for (std::size_t i=0; i<nsf+3*nvf; i++){
gfield0.at(info.current[0],i) += oldData[id][i] / info.previous.size();
gfieldAdapt.at(info.current[0],i) += oldData[id][i] / info.previous.size();
}
gfield0.at(info.current[0], nsf+3*nvf) = std::max(gfield0.at(info.current[0], nsf+3*nvf), oldData[id][nsf+3*nvf]);
gfieldAdapt.at(info.current[0], nsf+3*nvf) = std::max(gfieldAdapt.at(info.current[0], nsf+3*nvf), oldData[id][nsf+3*nvf]);
}
}
}
meshr->adaptionCleanup();
}
targetCells.clear();
for (Cell & cell : meshr->getCells()){
long id = cell.getId();
if (gfield0.at(id, nsf+3*nvf))
targetCells.insert(id);
}
pod.setStaticMesh(false);
pod.reconstructFields(gfield0, meshr, fields, &targetCells);
if (rank==0){
std::cout<< " " << std::endl;
std::cout<< ">> Reconstruction coeffs hybrid interface:" << std::endl;
for (std::size_t i = 0; i < nf; ++i)
std::cout<< pod.getReconstructionCoeffs()[i] << std::endl;
std::cout<< " " << std::endl;
std::cout<< ">> Reconstruction error hybrid interface" << std::endl;
double maxerr = 0.0;
for (Cell & cell : meshr->getCells()){
long id = cell.getId();
for (std::size_t i = 0; i < nsf+3*nvf; ++i)
maxerr = std::max(maxerr, std::abs(gfield0.at(id, i) - gfieldAdapt.at(id, i)));
}
std::cout<< ">> Max error :" << maxerr << std::endl;
}
}
int main(int argc, char *argv[])
{
#if BITPIT_ENABLE_MPI
MPI_Init(&argc,&argv);
#endif
int nProcs;
int rank;
#if BITPIT_ENABLE_MPI
MPI_Comm_size(MPI_COMM_WORLD, &nProcs);
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
#else
nProcs = 1;
rank = 0;
#endif
try {
run(rank, nProcs);
} catch (const std::exception &exception) {
log::cout() << exception.what();
exit(1);
}
#if BITPIT_ENABLE_MPI
MPI_Finalize();
#endif
}
The Cell class defines the cells.
Definition cell.hpp:42
Input binary archive.
std::istream & getStream()
The POD (Proper Orthogonal Decomposition) class provides an interface for defining POD object.
Definition pod.hpp:41
std::size_t getModeCount()
Definition pod.cpp:315
void run()
Definition pod.cpp:836
void setMeshType(MeshType type)
Definition pod.cpp:390
void setStaticMesh(bool flag)
Definition pod.cpp:483
void setSensorMask(const PiercedStorage< bool > &mask, VolumeKernel *mesh=nullptr)
Definition pod.cpp:672
void setMemoryMode(MemoryMode mode)
Definition pod.cpp:505
void setDirectory(const std::string &directory)
Definition pod.cpp:180
void reconstructFields(pod::PODField &field, pod::PODField &recon)
Definition pod.cpp:1631
void setWriteMode(WriteMode mode)
Definition pod.cpp:576
void restore()
Definition pod.cpp:864
void setModeCount(std::size_t nmodes)
Definition pod.cpp:304
void setName(const std::string &name)
Definition pod.cpp:160
std::vector< std::vector< double > > getReconstructionCoeffs()
Definition pod.cpp:800
void addSnapshot(const std::string &directory, const std::string &name)
Definition pod.cpp:208
void markCellForRefinement(long id)
PiercedVector< Cell > & getCells()
std::vector< long > findCellNeighs(long id) const
std::vector< adaption::Info > adaptionPrepare(bool trackAdaption=true)
void markCellForCoarsening(long id)
void restore(std::istream &stream, bool reregister=false)
virtual long getCellCount() const
std::vector< adaption::Info > adaptionAlter(bool trackAdaption=true, bool squeezeStorage=false)
Metafunction for generating a pierced storage.
__PS_REFERENCE__ at(id_t id, std::size_t k=0)
Metafunction for generating a pierced vector.
iterator insert(id_t id, const value_t &value)
The VolOctree defines a Octree patch.
Definition voloctree.hpp:37
The VolumeKernel class provides an interface for defining volume patches.
The Info struct defines the information associated to an adaption.
Definition adaption.hpp:63
The PODfield structure is used to store the fields inside POD classes.
std::unique_ptr< pod::ScalarStorage > scalar
std::unique_ptr< PiercedStorage< bool > > mask
std::unique_ptr< pod::VectorStorage > vector
--- layout: doxygen_footer ---