30#include <unordered_map>
31#include <unordered_set>
37#include "bitpit_private_lapacke.hpp"
38#include "bitpit_voloctree.hpp"
41#include "pod_voloctree.hpp"
71: m_filter(1), m_sensorMask(1)
75 m_communicator = MPI_COMM_NULL;
78 m_podkernel =
nullptr;
84 m_nModes = std::numeric_limits<std::size_t>::max();
88 m_nReconstructionSnapshots = 0;
94 initializeCommunicator(comm);
95 MPI_Comm_size(m_communicator, &m_nProcs);
96 MPI_Comm_rank(m_communicator, &m_rank);
132 m_correlationMatrices.clear();
134 m_minimizationMatrices.clear();
137 m_reconstructionDatabase.clear();
138 m_leave1outOffDatabase.clear();
144 m_nReconstructionSnapshots = 0;
147 m_podkernel->clearMapper();
150# if BITPIT_ENABLE_MPI
183 if (stat(directory.c_str(), &info) != 0) {
184 throw std::runtime_error(
"The directory \"" + directory +
"\" does not exists!");
185 }
else if (!(info.st_mode & S_IFDIR)) {
186 throw std::runtime_error(
"The directory \"" + directory +
"\" is not valid: a file with the same name was found!");
189 m_directory = directory +
"/";
221 m_database.push_back(file);
232 m_database = database;
233 m_nSnapshots = database.size();
256 m_leave1outOffDatabase.push_back(file);
264 std::size_t nsoff=m_leave1outOffDatabase.size();
265 for (std::size_t i=0; i<m_nSnapshots; ++i){
266 std::string snapi=m_database[i].directory+
"/"+m_database[i].name;
267 for (std::size_t j=0; j<nsoff; ++j){
268 std::string snapj=m_leave1outOffDatabase[j].directory+
"/"+m_leave1outOffDatabase[j].name;
269 if (m_database[i].directory == m_leave1outOffDatabase[j].directory && m_database[i].name == m_leave1outOffDatabase[j].name)
270 m_listActiveIDsLeave1out.push_back(i);
295 m_reconstructionDatabase.push_back(file);
296 m_nReconstructionSnapshots++;
328 m_energyLevel = energy - ENERGY_CHECK_TOLERANCE;
338 return m_energyLevel + ENERGY_CHECK_TOLERANCE;
348 m_errorThreshold = threshold;
358 return m_errorThreshold;
369 std::size_t nsf= namesf.size();
370 std::size_t nvf= namevf.size();
372 m_nameTargetErrorFields.clear();
374 for (std::size_t ifield=0; ifield<nsf; ++ifield)
375 m_nameTargetErrorFields[namesf[ifield]] = ifield;
377 for (std::size_t ifield=0; ifield<nvf; ++ifield){
378 for (std::size_t j=0; j<3; ++j)
379 m_nameTargetErrorFields[namevf[ifield][j]] = ifield*3+nsf+j;
392 throw std::runtime_error (
"POD mesh type already set. Change not allowed.");
400 m_podkernel = std::unique_ptr<PODKernel>(
new PODVolOctree(m_communicator));
402 m_podkernel = std::unique_ptr<PODKernel>(
new PODVolOctree());
407 throw std::runtime_error (
"POD mesh type not allowed.");
421 throw std::runtime_error (
"POD mesh type not set.");
424 if (_octreecast_mesh){
426 throw std::runtime_error (
"POD mesh type not set to VolOctree.");
427 m_podkernel->setMesh(std::move(mesh));
430 throw std::runtime_error (
"POD mesh type not allowed.");
442void POD::setMesh(
const std::string &directory,
const std::string &name)
458 throw std::runtime_error (
"POD mesh type not set.");
460 std::unique_ptr<VolumeKernel> mesh = m_podkernel->readMesh(file);
506 if (m_memoryMode == mode) {
543 if (m_runMode == mode) {
577 if (m_writeMode == mode)
599 if (m_reconstructionMode == mode)
624 return m_reconstructionMode;
636 if (m_errorMode == mode)
673 m_sensorMask.unsetKernel();
674 m_sensorMask.setStaticKernel(&(m_podkernel->getMesh()->getCells()));
676 if (m_staticMesh || mesh == m_podkernel->getMesh()){
678 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
679 long id = cell.getId();
680 m_sensorMask[id] = mask[id];
687 throw std::runtime_error (
"POD: null mesh pointer passed to setSensorMask");
690 if (m_podkernel->isMapperDirty())
691 _computeMapper(mesh);
693 std::unordered_set<long> trueCells;
695 long id = cell.getId();
697 trueCells.insert(
id);
701 std::unordered_set<long> mappedCells = m_podkernel->mapCellsToPOD(&trueCells);
703 m_podkernel->mapBoolFieldToPOD(mask, mesh, &mappedCells, m_sensorMask);
729 return m_nameScalarFields;
739 return m_nameVectorFields;
750 std::vector<std::string> names;
751 for (
const std::string &ss : m_nameScalarFields)
753 for (
const std::array<std::string,3> &ssa : m_nameVectorFields){
754 for (
const std::string &ss : ssa)
768 if (m_podkernel ==
nullptr)
769 throw std::runtime_error (
"POD mesh not built");
771 return m_podkernel->getMesh();
801 return m_reconstructionCoeffs;
809 return m_listActiveIDs;
817 return m_sizeInternal;
839 log::cout() <<
"pod : evaluate decomposition ... " << std::endl;
842 log::cout() <<
"pod : restore decomposition ... " << std::endl;
847 log::cout() <<
"pod : evaluate reconstruction ... " << std::endl;
875 log::cout() <<
"pod : computing mean field and pod mesh... " << std::endl;
882 log::cout() <<
"pod : computing correlation matrix... " << std::endl;
886 log::cout() <<
"pod : computing eigenvectors... " << std::endl;
890 log::cout() <<
"pod : computing pod modes... " << std::endl;
895 log::cout() <<
"pod : dumping... " << std::endl;
905 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
907 log::cout() <<
"Dumping POD ... " << std::endl;
909 std::string header =
"pod";
910 OBinaryArchive binaryWriter(m_directory + m_name, ARCHIVE_VERSION, header, dumpBlock);
911 std::ostream &dataStream = binaryWriter.
getStream();
915 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
921 binaryWriter.close();
925 for (std::size_t ir = 0; ir < m_nModes; ++ir)
937 std::string header =
"podmask";
938 OBinaryArchive binaryWriter(m_directory + m_name +
".mask.data", ARCHIVE_VERSION, header, dumpBlock);
939 std::ostream &dataStream = binaryWriter.
getStream();
940 m_sensorMask.dump(dataStream);
941 binaryWriter.close();
946 std::string header =
"podmesh";
947 OBinaryArchive binaryWriter(m_directory + m_name +
".mesh", ARCHIVE_VERSION, header, dumpBlock);
948 m_podkernel->getMesh()->dump(binaryWriter.
getStream());
949 binaryWriter.close();
954 std::vector<std::string > datastring;
955 m_podkernel->getMesh()->getVTK().setDirectory(m_directory);
956 m_podkernel->getMesh()->getVTK().setName(m_name);
958 std::size_t nf=m_nModes;
961 std::vector<std::vector<std::vector<double>>> fields(nf, std::vector<std::vector<double>> (m_nScalarFields, std::vector<double>(m_podkernel->getMesh()->getInternalCellCount(), 0.0)));
962 std::vector<std::vector<std::vector<std::array<double,3>>>> fieldv(nf, std::vector<std::vector<std::array<double,3>>>(m_nVectorFields, std::vector<std::array<double,3>>(m_podkernel->getMesh()->getInternalCellCount(), {{0.0, 0.0, 0.0}})));
965 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
967 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
968 if (cell.isInterior()) {
969 long id = cell.getId();
970 fields[0][isf][i] = m_mean.scalar->at(
id, isf);
974 m_podkernel->getMesh()->getVTK().addData(
"scalarField" + std::to_string(isf) +
"_mean", VTKFieldType::SCALAR, VTKLocation::CELL, fields[0][isf]);
975 datastring.push_back(
"scalarField" + std::to_string(isf) +
"_mean");
977 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
979 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
980 if (cell.isInterior()) {
981 long id = cell.getId();
982 fieldv[0][ivf][i] = m_mean.vector->at(
id, ivf);
986 m_podkernel->getMesh()->getVTK().addData(
"vectorField"+std::to_string(ivf)+
"_mean", VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[0][ivf]);
987 datastring.push_back(
"vectorField"+std::to_string(ivf)+
"_mean");
990 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
994 std::size_t nir = ir;
998 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
1000 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
1001 if (cell.isInterior()) {
1002 long id = cell.getId();
1003 fields[nir][isf][i] = m_modes[ir].scalar->at(
id, isf);
1007 m_podkernel->getMesh()->getVTK().addData(
"scalarField"+std::to_string(isf)+
"_podMode"+std::to_string(ir), VTKFieldType::SCALAR, VTKLocation::CELL, fields[nir][isf]);
1008 datastring.push_back(
"scalarField"+std::to_string(isf)+
"_podMode"+std::to_string(ir));
1011 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
1013 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
1014 if (cell.isInterior()) {
1015 long id = cell.getId();
1016 fieldv[nir][ivf][i] = m_modes[ir].vector->at(
id, ivf);
1020 m_podkernel->getMesh()->getVTK().addData(
"vectorField"+std::to_string(ivf)+
"_podMode"+std::to_string(ir), VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[nir][ivf]);
1021 datastring.push_back(
"vectorField"+std::to_string(ivf)+
"_podMode"+std::to_string(ir));
1025 m_modes[ir].clear();
1028 std::vector<uint8_t> mask(m_podkernel->getMesh()->getInternalCellCount(), 0);
1030 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
1031 if (cell.isInterior()) {
1032 long id = cell.getId();
1033 mask[i] = m_filter[id];
1038 m_podkernel->getMesh()->getVTK().addData(
"filter", VTKFieldType::SCALAR, VTKLocation::CELL, mask);
1039 datastring.push_back(
"filter");
1040 m_podkernel->getMesh()->write();
1042 for (std::size_t is = 0; is < datastring.size(); ++is)
1043 m_podkernel->getMesh()->getVTK().removeData(datastring[is]);
1052 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
1055 log::cout() <<
"Restore POD ... " << std::endl;
1057 std::string header =
"pod";
1059 std::istream &dataStream = binaryReader.
getStream();
1065 m_nModes = std::min(m_nModes, nr_);
1066 m_modes.resize(m_nModes);
1068 m_database.resize(m_nSnapshots);
1069 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
1075 binaryReader.close();
1081 podfile.
name = m_name;
1082 m_podkernel->restoreMesh(podfile);
1087 m_filter.unsetKernel();
1088 m_filter.setStaticKernel(&m_podkernel->getMesh()->getCells());
1089 m_sensorMask.unsetKernel();
1090 m_sensorMask.setStaticKernel(&m_podkernel->getMesh()->getCells());
1093 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
1106 IBinaryArchive binaryReader(m_directory + m_name +
".mask.data", dumpBlock);
1107 std::istream &dataStream = binaryReader.
getStream();
1108 m_sensorMask.restore(dataStream);
1109 binaryReader.close();
1127 log::cout() <<
"pod : computing pod mesh... " << std::endl;
1134 log::cout() <<
"pod : computing correlation matrix... " << std::endl;
1141 std::vector<std::vector<double>> h_correlationMatrices=m_correlationMatrices;
1142 std::vector<pod::SnapshotFile> h_database=m_database;
1143 std::size_t h_nSnapshots=m_nSnapshots;
1144 std::size_t nl1o=m_listActiveIDsLeave1out.size();
1146 m_nSnapshots=h_nSnapshots-1;
1148 for (std::size_t i=0; i<nl1o; ++i){
1149 std::size_t
id = m_listActiveIDsLeave1out[i];
1150 log::cout() <<
">> leave-1-out : removing snapshot " << m_database[id].directory+
"/"+m_database[id].name<< std::endl;
1152 m_reconstructionDatabase.clear();
1153 m_nReconstructionSnapshots = 0;
1155 m_database.erase(m_database.begin()+
id);
1157 for (std::size_t ifield = 0; ifield<m_nFields; ++ifield){
1159 m_correlationMatrices[ifield].erase(m_correlationMatrices[ifield].begin()+
id*h_nSnapshots, m_correlationMatrices[ifield].begin()+
id*h_nSnapshots+h_nSnapshots);
1163 for (std::size_t j=0; j<m_nSnapshots; ++j){
1164 m_correlationMatrices[ifield].erase(m_correlationMatrices[ifield].begin()+it);
1165 it= it-1+h_nSnapshots;
1170 log::cout() <<
"pod : computing eigenvectors... " << std::endl;
1174 log::cout() <<
"pod : computing pod modes... " << std::endl;
1181 m_correlationMatrices=h_correlationMatrices;
1182 m_database=h_database;
1188 std::string ename = m_name +
".error";
1191 m_nSnapshots=h_nSnapshots;
1201 if (m_database.empty())
1202 throw std::runtime_error (
"POD database empty");
1208#if BITPIT_ENABLE_MPI
1209 m_podkernel->getMesh()->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_INTERNAL);
1211 m_podkernel->getMesh()->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_ALL);
1218void POD::_evalMeanMesh()
1221 if (m_podkernel->getMesh() ==
nullptr){
1224 std::unique_ptr<VolumeKernel> readMesh = m_podkernel->readMesh(m_database[0]);
1225 m_podkernel->setMesh(std::move(readMesh));
1230 for (std::size_t i = 1; i < m_nSnapshots; ++i) {
1231 log::cout() <<
"pod : evaluation POD mesh - use snapshot " << i+1 <<
"/" << m_nSnapshots << std::endl;
1232 std::unique_ptr<VolumeKernel> readMesh = m_podkernel->readMesh(m_database[i]);
1233 m_podkernel->adaptMeshToMesh(m_podkernel->getMesh(), readMesh.get());
1234 m_podkernel->clearMapper();
1241 pod::PODField readf;
1242 readSnapshot(m_database[0], readf);
1246 m_podkernel->evalCellsVolume();
1248 m_filter.unsetKernel();
1249 m_filter.setStaticKernel(&(m_podkernel->getMesh()->getCells()));
1250 m_filter.fill(
true);
1254 m_mean.scalar = std::unique_ptr<pod::ScalarStorage>(
new pod::ScalarStorage(m_nScalarFields, &m_podkernel->getMesh()->getCells()));
1255 m_mean.vector = std::unique_ptr<pod::VectorStorage>(
new pod::VectorStorage(m_nVectorFields, &m_podkernel->getMesh()->getCells()));
1256 m_mean.scalar->fill(0.0);
1257 m_mean.vector->fill(std::array<double, 3>{{0.0, 0.0, 0.0}});
1260 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
1261 log::cout() <<
"pod : evaluation mean - use snapshot " << i+1 <<
"/" << m_nSnapshots << std::endl;
1262 pod::PODField readf;
1263 readSnapshot(m_database[i], readf);
1266 _computeMapper(readf.mesh);
1267 readf = pod::PODField(m_podkernel->mapPODFieldToPOD(readf,
nullptr));
1268 m_podkernel->clearMapper();
1271 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
1272 m_filter[cell.getId()] = (m_filter[cell.getId()]) && (readf.mask->at(cell.getId()));
1273 if (m_filter[cell.getId()]){
1274 for (std::size_t j = 0; j < m_nScalarFields; ++j)
1275 m_mean.scalar->at(cell.getId(),j) += readf.scalar->at(cell.getId(),j) / double(m_nSnapshots);
1276 for (std::size_t j = 0; j < m_nVectorFields; ++j)
1277 m_mean.vector->at(cell.getId(),j) += readf.vector->at(cell.getId(),j) / double(m_nSnapshots);
1280 for (std::size_t j = 0; j < m_nScalarFields; ++j)
1281 m_mean.scalar->at(cell.getId(),j) = 0.0;
1282 for (std::size_t j = 0; j < m_nVectorFields; ++j)
1283 m_mean.vector->at(cell.getId(),j) = {{0.0, 0.0, 0.0}};
1290 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
1291 pod::PODField readf;
1292 readSnapshot(m_database[i], readf);
1295 _computeMapper(readf.mesh);
1296 readf = pod::PODField(m_podkernel->mapPODFieldToPOD(readf,
nullptr));
1297 m_podkernel->clearMapper();
1300 for (
const Cell &cell : m_podkernel->getMesh()->getCells())
1301 m_filter[cell.getId()] = (m_filter[cell.getId()]) && (readf.mask->at(cell.getId()));
1306 setSensorMask(m_filter, m_podkernel->getMesh());
1319 m_listActiveIDs.clear();
1320 std::unordered_set<long> listGhost;
1321 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
1322 long id = cell.getId();
1325 if (mask[
id] && cell.isInterior()) {
1326 if (cell.isInterior()) {
1327 m_listActiveIDs.insert(
id);
1329 listGhost.insert(
id);
1334 m_sizeInternal = m_listActiveIDs.size();
1347 for (std::size_t i = 0; i < m_nSnapshots; ++i){
1351 _computeMapper(snapi.
mesh);
1352 snapi =
pod::PODField(m_podkernel->mapPODFieldToPOD(snapi,
nullptr));
1353 m_podkernel->clearMapper();
1357 diff(&snapi, m_mean);
1359 for (std::size_t j = i; j < m_nSnapshots; ++j){
1360 log::cout() <<
"pod : evaluation " << i <<
"," << j <<
" term of correlation matrix " << std::endl;
1364 _computeMapper(snapj.
mesh);
1365 snapj =
pod::PODField(m_podkernel->mapPODFieldToPOD(snapj,
nullptr));
1366 m_podkernel->clearMapper();
1370 diff(&snapj, m_mean);
1372 evalCorrelationTerm(i, snapi, j, snapj);
1377# if BITPIT_ENABLE_MPI
1378 for (std::size_t i = 0; i < m_nFields; i++ )
1379 MPI_Allreduce(MPI_IN_PLACE, m_correlationMatrices[i].data(), m_nSnapshots*m_nSnapshots, MPI_DOUBLE, MPI_SUM, m_communicator);
1387void POD::initCorrelation()
1389 m_correlationMatrices.clear();
1390 m_correlationMatrices.resize(m_nFields,std::vector<double>(m_nSnapshots*m_nSnapshots, 0.0));
1396void POD::initErrorMaps()
1399 m_errorMap.mesh = m_podkernel->getMesh();
1400 m_errorMap.scalar = std::unique_ptr<pod::ScalarStorage>(
new pod::ScalarStorage(m_nScalarFields, &m_podkernel->getMesh()->getCells()));
1401 m_errorMap.vector = std::unique_ptr<pod::VectorStorage>(
new pod::VectorStorage(m_nVectorFields, &m_podkernel->getMesh()->getCells()));
1402 m_errorMap.scalar->fill(0.0);
1403 m_errorMap.vector->fill(std::array<double, 3>{{0.0, 0.0, 0.0}});
1404 m_errorMap.mask = std::unique_ptr<PiercedStorage<bool>>(
new PiercedStorage<bool>(1, &m_podkernel->getMesh()->getCells()));
1405 m_errorMap.mask->fill(0);
1406 for (
long id : m_listActiveIDs)
1407 m_errorMap.mask->at(
id)=1;
1414void POD::evalCorrelationTerm(
int i, pod::PODField & snapi,
int j, pod::PODField & snapj)
1417 std::unordered_set<long>::iterator it;
1418 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it){
1420 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
1421 if (m_nScalarFields){
1422 double* datasi = snapi.scalar->rawData(rawIndex);
1423 double* datasj = snapj.scalar->rawData(rawIndex);
1424 for (std::size_t ifield = 0; ifield < m_nScalarFields; ++ifield){
1425 m_correlationMatrices[ifield][i*m_nSnapshots+j] += (*datasi)*(*datasj)*getRawCellVolume(rawIndex);
1430 if (m_nVectorFields){
1431 std::array<double,3>* datavi = snapi.vector->rawData(rawIndex);
1432 std::array<double,3>* datavj = snapj.vector->rawData(rawIndex);
1433 for (std::size_t ifield = m_nScalarFields; ifield < m_nFields; ++ifield){
1434 m_correlationMatrices[ifield][i*m_nSnapshots+j] +=
dotProduct((*datavi),(*datavj))*getRawCellVolume(rawIndex);
1450 m_reconstructionCoeffs.clear();
1451 m_reconstructionCoeffs.resize(m_nFields, std::vector<double>(m_nModes,0.0));
1453 for (std::size_t i = 0; i < m_nReconstructionSnapshots; ++i){
1460 _computeMapper(snapi.
mesh);
1463 snapi =
pod::PODField(m_podkernel->mapPODFieldToPOD(snapi,
nullptr));
1464 m_podkernel->clearMapper();
1469 std::string rname = m_name +
".recon." + m_reconstructionDatabase[i].name;
1476 std::ofstream outr(m_directory + m_name +
".recon.coeff."+ m_reconstructionDatabase[i].name +
".dat", std::ofstream::out);
1477 for (std::size_t j = 0; j < m_nFields; ++j)
1478 outr << m_reconstructionCoeffs[j] << std::endl;
1484 std::string ename = m_name +
".error."+ m_reconstructionDatabase[i].name;
1490 std::ofstream oute(m_directory + m_name +
".error.norm."+ m_reconstructionDatabase[i].name +
".dat", std::ofstream::out);
1491 for (std::size_t j = 0; j < m_nFields; ++j)
1492 oute << enorm[j] << std::endl;
1506 log::cout()<<
"Reading error field ..."<< std::endl;
1513 std::vector<std::size_t> scalarIds, vectorIds;
1514 std::map<std::string, std::size_t> targetErrorFields = m_nameTargetErrorFields;
1517 std::size_t count=0;
1518 for (
const std::string &val : m_nameScalarFields) {
1519 std::map<std::string, std::size_t>::iterator found = targetErrorFields.find(val);
1520 if (found != targetErrorFields.end()) {
1521 scalarIds.push_back(count);
1522 targetErrorFields.erase(found);
1528 for (
const std::array<std::string,3> &valv : m_nameVectorFields) {
1530 std::array<std::string,3> toerase;
1531 for (
const std::string &val : valv) {
1532 std::map<std::string, std::size_t>::iterator found = targetErrorFields.find(val);
1533 if (found != targetErrorFields.end()) {
1539 vectorIds.push_back(count);
1540 for (std::size_t i = 0; i < 3; ++i)
1541 targetErrorFields.erase(toerase[i]);
1547 if (targetErrorFields.size() != 0) {
1548 std::string verror =
" ";
1549 for (
const auto &targetErrorField : targetErrorFields)
1550 verror += targetErrorField.first +
" ";
1552 throw std::runtime_error (
"POD: fields not found in POD base or incomplete vector: " + verror);
1555 std::size_t nsf = scalarIds.size();
1556 std::size_t nvf = vectorIds.size();
1558 std::vector<std::array<double, 3>> minBoxes, maxBoxes;
1559 minBoxes.resize(m_nFields, {{0.0, 0.0, 0.0}});
1560 maxBoxes.resize(m_nFields, {{0.0, 0.0, 0.0}});
1562 std::unordered_set<long>::iterator it;
1563 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it) {
1565 std::array<double, 3> cellCentroid = error.mesh->evalCellCentroid(
id);
1566 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
1569 double* datas = error.scalar->rawData(rawIndex);
1570 for (std::size_t i = 0; i < nsf; ++i){
1571 double* datasi = datas+scalarIds[i];
1572 if (*datasi >= m_errorThreshold){
1573 maxBoxes[scalarIds[i]]=
max(cellCentroid,maxBoxes[scalarIds[i]]);
1574 minBoxes[scalarIds[i]]=
min(cellCentroid,minBoxes[scalarIds[i]]);
1579 std::array<double,3>* datav = error.vector->rawData(rawIndex);
1580 for (std::size_t i = 0; i< nvf; ++i) {
1581 std::array<double,3>* datavi = datav+vectorIds[i];
1582 if (std::sqrt(
dotProduct((*datavi),(*datavi)) ) >= m_errorThreshold){
1583 maxBoxes[vectorIds[i]]=
max(cellCentroid,maxBoxes[vectorIds[i]]);
1584 minBoxes[vectorIds[i]]=
min(cellCentroid,minBoxes[vectorIds[i]]);
1590# if BITPIT_ENABLE_MPI
1591 MPI_Allreduce(MPI_IN_PLACE,
reinterpret_cast<double *
>(maxBoxes.data()), m_nFields*3, MPI_DOUBLE, MPI_MAX, m_communicator);
1592 MPI_Allreduce(MPI_IN_PLACE,
reinterpret_cast<double *
>(minBoxes.data()), m_nFields*3, MPI_DOUBLE, MPI_MIN, m_communicator);
1595 log::cout()<<
">> Fields Bounding boxes " << std::endl;
1596 log::cout()<<
"min: "<< minBoxes << std::endl;
1597 log::cout()<<
"max: "<< maxBoxes << std::endl;
1599 std::array<double, 3> minBox = {{0.0, 0.0, 0.0}};
1600 std::array<double, 3> maxBox = {{0.0, 0.0, 0.0}};
1602 for (std::size_t i =0; i < m_nFields; ++i){
1603 maxBox=
max(maxBoxes[i],maxBox);
1604 minBox=
min(minBoxes[i],minBox);
1607 log::cout()<<
">> Bounding box " << std::endl;
1608 log::cout()<<
"("<< minBox <<
") ("<< maxBox <<
")"<< std::endl;
1611#if BITPIT_ENABLE_MPI
1612 runSolver = (m_rank == 0);
1617 std::ofstream outBB(m_directory + m_name +
".box.dat", std::ofstream::out);
1618 outBB << minBox << std::endl;
1619 outBB << maxBox << std::endl;
1636 evalReconstructionCoeffs(field);
1637 buildFields(m_reconstructionCoeffs, reconi);
1641 buildErrorMaps(field,reconi);
1654 reconi.
setMesh(m_podkernel->getMesh());
1655 buildFields(reconstructionCoeffs, reconi);
1667 std::map<std::string, std::size_t> targetFields,
1668 const std::unordered_set<long> * targetCells)
1670 std::vector<std::size_t> podscalarIds;
1671 std::vector<std::size_t> podvectorIds;
1673 std::vector<std::size_t> scalarIds;
1674 std::vector<std::array<std::size_t, 3>> vectorIds;
1677 std::size_t count = 0;
1678 for (
const std::string &val : m_nameScalarFields) {
1679 std::map<std::string, std::size_t>::iterator found = targetFields.find(val);
1680 if (found != targetFields.end()) {
1681 std::size_t ind = found->second;
1682 podscalarIds.push_back(count);
1683 scalarIds.push_back(ind);
1684 targetFields.erase(found);
1691 for (
const std::array<std::string,3> &valv : m_nameVectorFields) {
1693 std::array<std::size_t,3> vectorarr;
1694 std::array<std::string,3> toerase;
1695 for (
const std::string &val : valv) {
1696 std::map<std::string, std::size_t>::iterator found = targetFields.find(val);
1697 if (found != targetFields.end()) {
1698 std::size_t ind = found->second;
1699 vectorarr[ic] = ind;
1705 podvectorIds.push_back(count);
1706 vectorIds.push_back(vectorarr);
1707 for (std::size_t i = 0; i < 3; ++i)
1708 targetFields.erase(toerase[i]);
1714 if (targetFields.size() != 0) {
1715 std::string verror =
" ";
1716 for (
const auto &targetField : targetFields)
1717 verror += targetField.first +
" ";
1719 throw std::runtime_error (
"POD: fields not found in POD base or incomplete vector: " + verror);
1725 evalReconstructionCoeffs(fields, scalarIds, podscalarIds, vectorIds, podvectorIds);
1727 buildFields(m_reconstructionCoeffs, fields, scalarIds, podscalarIds, vectorIds, podvectorIds, targetCells);
1733 if (m_podkernel->isMapperDirty())
1734 _computeMapper(mesh);
1737 PiercedStorage<double> mappedFields = m_podkernel->mapFieldsToPOD(fields, mesh, &m_listActiveIDs, scalarIds, vectorIds);
1739 evalReconstructionCoeffs(mappedFields, scalarIds, podscalarIds, vectorIds, podvectorIds);
1742 std::unordered_set<long> mappedCells;
1743 std::unordered_set<long>* ptr_mappedCells =
nullptr;
1745 mappedCells = m_podkernel->mapCellsToPOD(targetCells);
1746 ptr_mappedCells = &mappedCells;
1749 buildFields(m_reconstructionCoeffs, mappedFields, scalarIds, podscalarIds, vectorIds, podvectorIds, ptr_mappedCells);
1751 m_podkernel->mapFieldsFromPOD(fields, mesh, targetCells, mappedFields, scalarIds, vectorIds);
1765 std::vector<std::vector<double>> reconstructionCoeffs(m_nFields, std::vector<double>(m_nModes,0.0));
1767 diff(&field, m_mean);
1771 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
1777 for (
long id : m_listActiveIDs) {
1778 std::size_t rawIndex = field.
mesh->
getCells().getRawIndex(
id);
1781 if (m_nScalarFields) {
1782 const double *modes = m_modes[ir].scalar->rawData(rawIndex);
1783 const double *data = field.
scalar->rawData(rawIndex);
1784 for (std::size_t ifs = 0; ifs < m_nScalarFields; ++ifs) {
1785 reconstructionCoeffs[ifs][ir] += data[ifs] * modes[ifs] * getRawCellVolume(rawIndex);
1790 if (m_nVectorFields) {
1791 const std::array<double,3> *modes = m_modes[ir].vector->rawData(rawIndex);
1792 const std::array<double,3> *data = field.
vector->rawData(rawIndex);
1793 for (std::size_t ifv = 0; ifv < m_nVectorFields; ++ifv) {
1794 reconstructionCoeffs[ifv+m_nScalarFields][ir] +=
dotProduct(data[ifv], modes[ifv]) * getRawCellVolume(rawIndex);
1800#if BITPIT_ENABLE_MPI
1801 for (std::size_t i = 0; i < m_nFields; ++i) {
1802 MPI_Allreduce(MPI_IN_PLACE, reconstructionCoeffs[i].data(), m_nModes, MPI_DOUBLE, MPI_SUM, m_communicator);
1807 sum(&field, m_mean);
1810 return reconstructionCoeffs;
1821 std::vector<double>
norm = fieldsMax(snap);
1823 for (
long id : m_listActiveIDs){
1824 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
1825 if (m_nScalarFields){
1826 double* recons = recon.
scalar->rawData(rawIndex);
1827 double* snaps = snap.
scalar->rawData(rawIndex);
1828 double* errors = m_errorMap.scalar->rawData(rawIndex);
1829 for (std::size_t ifs = 0; ifs < m_nScalarFields; ++ifs){
1830 *errors=std::max(*errors, std::abs(*recons - *snaps)/
norm[ifs]);
1836 if (m_nVectorFields){
1837 std::array<double,3>* reconv = recon.
vector->rawData(rawIndex);
1838 std::array<double,3>* snapv = snap.
vector->rawData(rawIndex);
1839 std::array<double,3>* errorv = m_errorMap.vector->rawData(rawIndex);
1840 for (std::size_t ifv = 0; ifv < m_nVectorFields; ++ifv){
1841 for (std::size_t j=0; j<3; ++j)
1842 (*errorv)[j]=std::max((*errorv)[j], std::abs((*reconv)[j] - (*snapv)[j])/
norm[ifv]);
1855void POD::initMinimization()
1857 m_minimizationMatrices.clear();
1858 m_minimizationMatrices.resize(m_nFields,std::vector<double>(m_nModes*m_nModes, 0.0));
1864void POD::evalMinimizationMatrices()
1870 if (m_reconstructionMode == ReconstructionMode::PROJECTION){
1871 for (std::size_t ifield = 0; ifield < m_nFields; ++ifield){
1872 for (std::size_t ir = 0; ir < m_nModes; ++ir)
1873 m_minimizationMatrices[ifield][ir*m_nModes+ir] = 1.0;
1878 for (std::size_t ir = 0; ir < m_nModes; ++ir){
1879 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
1882 for (std::size_t jr = ir; jr < m_nModes; ++jr){
1883 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
1886 std::unordered_set<long>::iterator it;
1887 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it){
1889 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
1890 if (m_nScalarFields){
1891 double* datasi = m_modes[ir].scalar->rawData(rawIndex);
1892 double* datasj = m_modes[jr].scalar->rawData(rawIndex);
1893 for (std::size_t ifield = 0; ifield < m_nScalarFields; ++ifield){
1894 m_minimizationMatrices[ifield][ir*m_nModes+jr] += (*datasi)*(*datasj)*getRawCellVolume(rawIndex);
1900 for (std::size_t ifield = 0; ifield < m_nScalarFields; ++ifield){
1901 m_minimizationMatrices[ifield][jr*m_nModes+ir] = m_minimizationMatrices[ifield][ir*m_nModes+jr];
1905 if (m_nVectorFields){
1906 std::array<double,3>* datavi = m_modes[ir].vector->rawData(rawIndex);
1907 std::array<double,3>* datavj = m_modes[jr].vector->rawData(rawIndex);
1908 for (std::size_t ifield = m_nScalarFields; ifield < m_nFields; ++ifield){
1909 m_minimizationMatrices[ifield][ir*m_nModes+jr] +=
dotProduct((*datavi),(*datavj))*getRawCellVolume(rawIndex);
1915 for (std::size_t ifield = m_nScalarFields; ifield < m_nFields; ++ifield)
1916 m_minimizationMatrices[ifield][jr*m_nModes+ir] = m_minimizationMatrices[ifield][ir*m_nModes+jr];
1923# if BITPIT_ENABLE_MPI
1924 for (std::size_t i = 0; i < m_nFields; i++ )
1925 MPI_Allreduce(MPI_IN_PLACE, m_minimizationMatrices[i].data(), m_nModes*m_nModes, MPI_DOUBLE, MPI_SUM, m_communicator);
1938void POD::solveMinimization(std::vector<std::vector<double> > & rhs)
1941#if BITPIT_ENABLE_MPI
1942 runSolver = (m_rank == 0);
1948 for (std::size_t i = 0; i < m_nFields; ++i) {
1949 std::vector<double> A(m_nModes * m_nModes);
1950 std::vector<int> ipiv(m_nModes);
1952 for (std::size_t j = 0; j < m_nModes*m_nModes; ++j)
1953 A[j] = m_minimizationMatrices[i][j];
1955 lapack_int
info = LAPACKE_dgesv(LAPACK_COL_MAJOR, m_nModes, 1, A.data(), m_nModes, ipiv.data(), rhs[i].data(), m_nModes);
1957 log::cout() <<
"WARNING! algorithm convergence info " <<
info << std::endl;
1958 throw std::runtime_error(
"Unable to solve minimization problem.");
1963 for (std::size_t i = 0; i < m_nFields; ++i) {
1964 rhs[i] = std::vector<double>(m_nModes, 0.0);
1975 m_lambda.resize(m_nFields,std::vector<double>(m_nSnapshots, 0.0));
1977 m_podCoeffs.clear();
1978 m_podCoeffs.resize(m_nFields, std::vector<std::vector<double>>(m_nSnapshots, std::vector<double>(m_nSnapshots,0.0)));
1981#if BITPIT_ENABLE_MPI
1982 runSolver = (m_rank == 0);
1988 int N = m_nSnapshots*m_nSnapshots;
1989 for (std::size_t i = 0; i < m_nFields; ++i) {
1990 std::vector<double> Marr(N);
1991 std::vector<double> alambda(m_nSnapshots);
1992 for (std::size_t j = 0; j < m_nSnapshots*m_nSnapshots; ++j)
1993 Marr[j/m_nSnapshots+j%m_nSnapshots*m_nSnapshots] = m_correlationMatrices[i][j];
1995 lapack_int info = LAPACKE_dsyev(LAPACK_COL_MAJOR,
'V',
'U', m_nSnapshots, Marr.data(), m_nSnapshots, alambda.data());
1997 log::cout() <<
"WARNING! algorithm convergence info " << info << std::endl;
1998 throw std::runtime_error(
"Unable to solve minimization problem.");
2002 checkModeCount(alambda.data(), i);
2004 for (std::size_t n = 0; n < m_nSnapshots; ++n)
2005 m_lambda[i][n] = alambda[m_nSnapshots - n - 1];
2007 for (std::size_t n = 0; n < m_nSnapshots; ++n) {
2008 for (std::size_t p = 0; p < m_nSnapshots; ++p)
2009 m_podCoeffs[i][p][n] = Marr[N - m_nSnapshots * (p + 1) + n] * std::sqrt(std::abs(m_lambda[i][p]));
2014 for (std::size_t ir=0; ir<m_nModes; ++ir) {
2015 std::ofstream out(m_directory + m_name +
".coeff."+std::to_string(ir)+
".dat", std::ofstream::out);
2016 for (std::size_t i = 0; i < m_nFields; ++i)
2017 out << m_podCoeffs[i][ir] << std::endl;
2022 std::ofstream outl(m_directory + m_name +
".lambda.dat", std::ofstream::out);
2023 for (std::size_t i = 0; i < m_nFields; ++i)
2024 outl << m_lambda[i] << std::endl;
2030#if BITPIT_ENABLE_MPI
2031 long bufferSize = m_nModes;
2032 MPI_Bcast(&bufferSize, 1, MPI_LONG, 0, m_communicator);
2033 m_nModes = bufferSize;
2034 for (std::size_t i = 0; i < m_nFields; ++i) {
2035 for (std::size_t p=0; p<m_nModes; ++p)
2036 MPI_Bcast(m_podCoeffs[i][p].data(), m_nSnapshots, MPI_DOUBLE, 0, m_communicator);
2038 MPI_Bcast(m_lambda[i].data(), m_nSnapshots, MPI_DOUBLE, 0, m_communicator);
2049void POD::checkModeCount(
double *alambda, std::size_t ifield)
2054 m_nModes = std::min(m_nModes, m_nSnapshots-1);
2056 m_nModes = std::min(m_nModes, m_nSnapshots);
2065 std::vector<double> dlambda(m_nSnapshots);
2067 for (std::size_t n=hn; n<m_nSnapshots; ++n){
2068 dlambda[n] = alambda[m_nSnapshots-n-1];
2069 sumen += std::abs(dlambda[n]);
2071 while (nr < m_nModes && en < m_energyLevel){
2072 en += std::abs(dlambda[nr]) / sumen * 100;
2075 _m_nr.push_back(nr);
2076 if (ifield == m_nFields-1)
2077 m_nModes = *(std::max_element(_m_nr.begin(), _m_nr.end()));
2085 m_modes.resize(m_nModes);
2092void POD::_evalModes()
2094 for (std::size_t ir = 0; ir < m_nModes; ++ir){
2095 m_modes[ir].scalar = std::unique_ptr<pod::ScalarStorage>(
new pod::ScalarStorage(m_nScalarFields, &m_podkernel->getMesh()->getCells()));
2096 m_modes[ir].vector = std::unique_ptr<pod::VectorStorage>(
new pod::VectorStorage(m_nVectorFields, &m_podkernel->getMesh()->getCells()));
2097 m_modes[ir].scalar->fill(0.0);
2098 m_modes[ir].vector->fill(std::array<double, 3>{{0.0, 0.0, 0.0}});
2101 for (std::size_t is = 0; is < m_nSnapshots; ++is){
2102 pod::PODField snapi;
2103 readSnapshot(m_database[is], snapi);
2105 _computeMapper(snapi.mesh);
2106 snapi = pod::PODField(m_podkernel->mapPODFieldToPOD(snapi,
nullptr));
2107 m_podkernel->clearMapper();
2111 diff(&snapi, m_mean);
2113 for (std::size_t ir = 0; ir < m_nModes; ++ir){
2114 for (
long id : m_listActiveIDs){
2115 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2116 if (m_nScalarFields){
2117 double* modes = m_modes[ir].scalar->rawData(rawIndex);
2118 double* datas = snapi.scalar->rawData(rawIndex);
2119 for (std::size_t ifs = 0; ifs < m_nScalarFields; ++ifs){
2120 *modes = *modes + *datas*m_podCoeffs[ifs][ir][is] / m_lambda[ifs][ir];
2125 if (m_nVectorFields){
2126 std::array<double, 3>* modev = m_modes[ir].vector->rawData(rawIndex);
2127 std::array<double, 3>* datav = snapi.vector->rawData(rawIndex);
2128 for (std::size_t ifv = 0; ifv < m_nVectorFields; ++ifv){
2129 *modev = *modev + *datav*m_podCoeffs[m_nScalarFields+ifv][ir][is] / m_lambda[m_nScalarFields+ifv][ir];
2135 if (m_memoryMode == MemoryMode::MEMORY_LIGHT){
2137 m_modes[ir].clear();
2152 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2154 std::string filename = std::string(snap.
directory) +
"/" + std::string(snap.
name) +
".data";
2156 std::istream &dataStream = binaryReader.
getStream();
2160 fieldr.
setMesh(m_podkernel->getMesh());
2164 fieldr.
setMesh(m_podkernel->readMesh(snap));
2170 fieldr.
mask->restore(dataStream);
2176 m_nameScalarFields.resize(m_nScalarFields);
2177 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2180 if (m_nScalarFields) {
2182 fieldr.
scalar->restore(dataStream);
2187 m_nFields = m_nScalarFields + m_nVectorFields;
2188 m_nameVectorFields.resize(m_nVectorFields);
2189 for (std::size_t i = 0; i < m_nVectorFields; ++i)
2192 if (m_nVectorFields) {
2194 fieldr.
vector->restore(dataStream);
2196 binaryReader.close();
2200 std::vector<std::string > datastring;
2204 std::vector<std::vector<double>> fields(m_nScalarFields, std::vector<double>(mesh->
getCellCount()));
2205 std::vector<std::vector<std::array<double,3>>> fieldv(m_nVectorFields, std::vector<std::array<double,3>>(mesh->
getCellCount()));
2209 long id = cell.getId();
2210 mask[i] = uint8_t(fieldr.
mask->at(
id));
2213 mesh->
getVTK().
addData(
"mask", VTKFieldType::SCALAR, VTKLocation::CELL, mask);
2214 datastring.push_back(
"mask");
2216 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
2219 long id = cell.getId();
2220 fields[isf][i] = fieldr.
scalar->at(
id, isf);
2223 mesh->
getVTK().
addData(m_nameScalarFields[isf], VTKFieldType::SCALAR, VTKLocation::CELL, fields[isf]);
2224 datastring.push_back(m_nameScalarFields[isf]);
2227 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
2230 long id = cell.getId();
2231 fieldv[ivf][i] = fieldr.
vector->at(
id, ivf);
2234 std::string vname= m_nameVectorFields[ivf][0].substr(0,m_nameVectorFields[ivf][0].size()-2);
2235 mesh->
getVTK().
addData(vname, VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[ivf]);
2236 datastring.push_back(vname);
2240 if (mesh == m_podkernel->getMesh())
2245 for (std::size_t is=0; is < datastring.size(); ++is)
2250 if (mesh == m_podkernel->getMesh()){
2251#if BITPIT_ENABLE_MPI
2264void POD::readMode(std::size_t ir)
2266 if (!m_podkernel->getMesh())
2267 throw std::runtime_error (
"POD: POD mesh not initialized before read a POD mode");
2270 throw std::runtime_error (
"POD: mode index > number of POD modes during read");
2272 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2274 if (ir == m_nModes){
2276 std::string filename = m_directory + m_name +
".mean"+
".data";
2278 std::istream &dataStream = binaryReader.getStream();
2281 m_filter.restore(dataStream);
2287 m_nameScalarFields.resize(m_nScalarFields);
2288 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2291 if (m_nScalarFields) {
2292 m_mean.scalar = std::unique_ptr<pod::ScalarStorage>(
new pod::ScalarStorage(m_nScalarFields, &m_podkernel->getMesh()->getCells()));
2293 m_mean.scalar->restore(dataStream);
2297 utils::binary::read(dataStream, m_nVectorFields);
2298 m_nameVectorFields.resize(m_nVectorFields);
2299 m_nFields = m_nScalarFields + m_nVectorFields;
2300 for (std::size_t i = 0; i < m_nVectorFields; ++i)
2301 utils::binary::read(dataStream, m_nameVectorFields[i]);
2303 if (m_nVectorFields) {
2304 m_mean.vector = std::unique_ptr<pod::VectorStorage>(
new pod::VectorStorage(m_nVectorFields, &m_podkernel->getMesh()->getCells()));
2305 m_mean.vector->restore(dataStream);
2307 binaryReader.close();
2311 m_modes[ir].clear();
2313 std::string filename = m_directory + m_name +
".mode."+std::to_string(ir)+
".data";
2314 IBinaryArchive binaryReader(filename, dumpBlock);
2315 std::istream &dataStream = binaryReader.getStream();
2318 m_filter.restore(dataStream);
2323 utils::binary::read(dataStream, m_nScalarFields);
2324 m_nameScalarFields.resize(m_nScalarFields);
2325 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2326 utils::binary::read(dataStream, m_nameScalarFields[i]);
2328 if (m_nScalarFields) {
2329 m_modes[ir].scalar = std::unique_ptr<pod::ScalarStorage>(
new pod::ScalarStorage(m_nScalarFields, &m_podkernel->getMesh()->getCells()));
2330 m_modes[ir].scalar->restore(dataStream);
2334 utils::binary::read(dataStream, m_nVectorFields);
2335 m_nFields = m_nScalarFields + m_nVectorFields;
2336 m_nameVectorFields.resize(m_nVectorFields);
2337 for (std::size_t i = 0; i < m_nVectorFields; ++i)
2338 utils::binary::read(dataStream, m_nameVectorFields[i]);
2340 if (m_nVectorFields) {
2341 m_modes[ir].vector = std::unique_ptr<pod::VectorStorage>(
new pod::VectorStorage(m_nVectorFields, &m_podkernel->getMesh()->getCells()));
2342 m_modes[ir].vector->restore(dataStream);
2344 binaryReader.close();
2353double POD::getCellVolume(
long id)
2355 return m_podkernel->getCellVolume(
id);
2363double POD::getRawCellVolume(
long rawIndex)
2365 return m_podkernel->getRawCellVolume(rawIndex);
2373void POD::dumpMode(std::size_t ir)
2376 throw std::runtime_error (
"POD: mode index > number of pod modes during write");
2378 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2380 if (ir == m_nModes) {
2382 std::string header =
"podmean." + std::to_string(ir);
2383 OBinaryArchive binaryWriter(m_directory + m_name +
".mean"+
".data", ARCHIVE_VERSION, header, dumpBlock);
2384 std::ostream &dataStream = binaryWriter.getStream();
2385 m_filter.dump(dataStream);
2386 utils::binary::write(dataStream,std::size_t(m_nScalarFields));
2387 if (m_nScalarFields) {
2388 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2389 utils::binary::write(dataStream, m_nameScalarFields[i]);
2390 m_mean.scalar->dump(dataStream);
2393 utils::binary::write(dataStream,std::size_t(m_nVectorFields));
2394 if (m_nVectorFields) {
2395 for (std::size_t i = 0; i < m_nVectorFields; ++i)
2396 utils::binary::write(dataStream, m_nameVectorFields[i]);
2397 m_mean.vector->dump(dataStream);
2399 binaryWriter.close();
2403 std::string header =
"podmode." + std::to_string(ir);
2404 OBinaryArchive binaryWriter(m_directory + m_name +
".mode."+std::to_string(ir)+
".data", ARCHIVE_VERSION, header, dumpBlock);
2405 std::ostream &dataStream = binaryWriter.getStream();
2406 m_filter.dump(dataStream);
2407 utils::binary::write(dataStream,std::size_t(m_nScalarFields));
2408 if (m_nScalarFields) {
2409 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2410 utils::binary::write(dataStream, m_nameScalarFields[i]);
2412 m_modes[ir].scalar->dump(dataStream);
2415 utils::binary::write(dataStream,std::size_t(m_nVectorFields));
2416 if (m_nVectorFields) {
2417 for (std::size_t i = 0; i < m_nVectorFields; ++i)
2418 utils::binary::write(dataStream, m_nameVectorFields[i]);
2420 m_modes[ir].vector->dump(dataStream);
2422 binaryWriter.close();
2434 log::cout() <<
"Dumping snapshot ... " << std::endl;
2436 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2440 std::string header = name;
2441 OBinaryArchive binaryWriter(m_directory + name +
".data", ARCHIVE_VERSION, header, dumpBlock);
2442 std::ostream &dataStream = binaryWriter.
getStream();
2443 m_filter.dump(dataStream);
2445 if (m_nScalarFields) {
2446 for (std::size_t i = 0; i < m_nScalarFields; ++i) {
2449 field.
scalar->dump(dataStream);
2453 if (m_nVectorFields) {
2454 for (std::size_t i = 0; i < m_nVectorFields; ++i) {
2457 field.
vector->dump(dataStream);
2459 binaryWriter.close();
2464 std::string header = name;
2465 OBinaryArchive binaryWriter(m_directory + name +
".mesh", ARCHIVE_VERSION, header, dumpBlock);
2466 std::ostream &dataStream = binaryWriter.
getStream();
2468 binaryWriter.close();
2486 std::size_t nsf = 0;
2487 if (m_nScalarFields){
2488 nsf = a.
scalar->getFieldCount();
2489 if (nsf != b.
scalar->getFieldCount())
2490 throw std::runtime_error (
"POD: different field count");
2493 std::size_t nvf = 0;
2494 if (m_nVectorFields){
2495 nvf = a.
vector->getFieldCount();
2496 if (nvf != b.
vector->getFieldCount())
2497 throw std::runtime_error (
"POD: different field count");
2501 for (
long id : a.
scalar->getKernel()->getIds()){
2502 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2503 double* datasa = a.
scalar->rawData(rawIndex);
2504 double* datasb = b.
scalar->rawData(rawIndex);
2505 for (std::size_t i = 0; i < nsf; ++i){
2506 *datasa = *datasa - *datasb;
2513 for (
long id : a.
scalar->getKernel()->getIds()){
2514 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2515 std::array<double,3>* datava = a.
vector->rawData(rawIndex);
2516 std::array<double,3>* datavb = b.
vector->rawData(rawIndex);
2517 for (std::size_t i = 0; i < nvf; ++i)
2518 *datava = *datava - *datavb;
2531void POD::sum(pod::PODField* _a,
const pod::PODMode& b)
2533 pod::PODField& a = *_a;
2534 std::size_t nsf = 0;
2535 if (m_nScalarFields){
2536 nsf = a.scalar->getFieldCount();
2537 if (nsf != b.scalar->getFieldCount())
2538 throw std::runtime_error (
"POD: different field count");
2541 std::size_t nvf = 0;
2542 if (m_nVectorFields){
2543 nvf = a.vector->getFieldCount();
2544 if (nvf != b.vector->getFieldCount())
2545 throw std::runtime_error (
"POD: different field count");
2549 for (
long id : a.scalar->getKernel()->getIds()){
2550 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2551 double *dataa = a.scalar->rawData(rawIndex);
2552 const double *datab = b.scalar->rawData(rawIndex);
2553 for (std::size_t i = 0; i < nsf; ++i){
2554 dataa[i] += datab[i];
2559 for (
long id : a.scalar->getKernel()->getIds()){
2560 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2561 std::array<double,3> *dataa = a.vector->rawData(rawIndex);
2562 const std::array<double,3> *datab = b.vector->rawData(rawIndex);
2563 for (std::size_t i = 0; i < nvf; ++i) {
2564 dataa[i] += datab[i];
2578 std::vector<double>
norm;
2579 norm.resize(m_nFields,0.0);
2581 std::unordered_set<long>::iterator it;
2582 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it){
2584 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2585 if (m_nScalarFields){
2586 double* datas = snap.
scalar->rawData(rawIndex);
2587 for (std::size_t ifield = 0; ifield < m_nScalarFields; ++ifield){
2588 norm[ifield] += (*datas)*(*datas)*snap.
mesh->evalCellVolume(
id);
2592 if (m_nVectorFields){
2593 std::array<double,3>* datav = snap.
vector->rawData(rawIndex);
2594 for (std::size_t ifield = m_nScalarFields; ifield < m_nFields; ++ifield){
2601# if BITPIT_ENABLE_MPI
2602 MPI_Allreduce(MPI_IN_PLACE,
norm.data(), m_nFields, MPI_DOUBLE, MPI_SUM, m_communicator);
2605 for (std::size_t ifield = 0; ifield < m_nFields; ++ifield)
2606 norm[ifield]=std::sqrt(
norm[ifield]);
2619 std::vector<double>
max;
2620 max.resize(m_nFields,0.0);
2621 for (
long id : m_listActiveIDs) {
2622 std::size_t rawIndex = snap.
mesh->
getCells().getRawIndex(
id);
2624 if (m_nScalarFields){
2625 const double *data = snap.
scalar->rawData(rawIndex);
2626 for (std::size_t i = 0; i < m_nScalarFields; ++i){
2627 double fieldData = data[i];
2628 max[i] = std::max(
max[i], std::abs(fieldData));
2631 if (m_nVectorFields){
2632 const std::array<double,3> *data = snap.
vector->rawData(rawIndex);
2633 for (std::size_t i = m_nScalarFields; i < m_nFields; ++i){
2634 const std::array<double,3> &fieldData = data[i - m_nScalarFields];
2640# if BITPIT_ENABLE_MPI
2641 MPI_Allreduce(MPI_IN_PLACE,
max.data(), m_nFields, MPI_DOUBLE, MPI_MAX, m_communicator);
2647#if BITPIT_ENABLE_MPI
2653void POD::initializeCommunicator(MPI_Comm communicator)
2656 if (isCommunicatorSet())
2657 throw std::runtime_error (
"POD communicator can be set just once");
2660 if (communicator == MPI_COMM_NULL)
2661 throw std::runtime_error (
"POD communicator is not valid");
2668 MPI_Comm_dup(communicator, &m_communicator);
2675MPI_Comm POD::getCommunicator()
const
2677 return m_communicator;
2686bool POD::isCommunicatorSet()
const
2689 return (getCommunicator() != MPI_COMM_NULL);
2695void POD::freeCommunicator()
2697 if (!isCommunicatorSet())
2700 int finalizedCalled;
2701 MPI_Finalized(&finalizedCalled);
2702 if (finalizedCalled)
2705 MPI_Comm_free(&m_communicator);
2714void POD::evalReconstructionCoeffs(pod::PODField & field)
2716 evalMinimizationMatrices();
2719 m_reconstructionCoeffs = projectField(field);
2722 solveMinimization(m_reconstructionCoeffs);
2734void POD::evalReconstructionCoeffs(PiercedStorage<double> &fields,
2735 const std::vector<std::size_t> &scalarIds,
const std::vector<std::size_t> & podscalarIds,
2736 const std::vector<std::array<std::size_t, 3>> & vectorIds,
const std::vector<std::size_t> & podvectorIds)
2738 std::size_t nsf = scalarIds.size();
2739 std::size_t nvf = vectorIds.size();
2742 diff(fields, m_mean, scalarIds, podscalarIds, vectorIds, podvectorIds, &m_listActiveIDs);
2744 m_reconstructionCoeffs.clear();
2745 m_reconstructionCoeffs.resize(m_nFields, std::vector<double>(m_nModes, 0.0));
2748 evalMinimizationMatrices();
2751 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
2752 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
2755 std::unordered_set<long>::iterator it;
2756 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it) {
2758 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2759 double *datag = fields.rawData(rawIndex);
2761 double *modes = m_modes[ir].scalar->rawData(rawIndex);
2762 for (std::size_t ifield = 0; ifield < nsf; ++ifield) {
2763 double *modesi = modes + podscalarIds[ifield];
2764 double *datagi = datag + scalarIds[ifield];
2765 m_reconstructionCoeffs[ifield][ir] += (*datagi)*(*modesi)*getRawCellVolume(rawIndex);
2769 std::array<double,3>* modev = m_modes[ir].vector->rawData(rawIndex);
2770 for (std::size_t ifield = 0; ifield < nvf; ++ifield) {
2771 std::array<double,3>* modevi = modev + podvectorIds[ifield];
2772 for (std::size_t j = 0; j < 3; ++j) {
2773 double *datagi = datag + vectorIds[ifield][j];
2774 m_reconstructionCoeffs[m_nScalarFields+ifield][ir] += ((*datagi)*(*modevi)[j])*getRawCellVolume(rawIndex);
2781#if BITPIT_ENABLE_MPI
2782 for (std::size_t i = 0; i < m_nFields; ++i)
2783 MPI_Allreduce(MPI_IN_PLACE, m_reconstructionCoeffs[i].data(), m_nModes, MPI_DOUBLE, MPI_SUM, m_communicator);
2787 solveMinimization(m_reconstructionCoeffs);
2790 sum(fields, m_mean, scalarIds, podscalarIds, vectorIds, podvectorIds, &m_listActiveIDs);
2800void POD::buildFields(
const std::vector<std::vector<double>> &reconstructionCoeffs, pod::PODField &recon)
2803 recon.scalar = std::unique_ptr<pod::ScalarStorage>(
new pod::ScalarStorage(m_nScalarFields, &m_podkernel->getMesh()->getCells()));
2804 recon.vector = std::unique_ptr<pod::VectorStorage>(
new pod::VectorStorage(m_nVectorFields, &m_podkernel->getMesh()->getCells()));
2805 recon.scalar->fill(0.0);
2806 recon.vector->fill(std::array<double, 3>{{0.0, 0.0, 0.0}});
2809 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
2810 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2815 for (
long id : m_listActiveIDs) {
2816 std::size_t rawIndex = recon.mesh->getCells().getRawIndex(
id);
2819 if (m_nScalarFields) {
2820 const double *modes = m_modes[ir].scalar->rawData(rawIndex);
2821 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
2822 recon.scalar->at(
id, isf) += modes[isf] * reconstructionCoeffs[isf][ir];
2827 if (m_nVectorFields) {
2828 const std::array<double,3> *modes = m_modes[ir].vector->rawData(rawIndex);
2829 for (std::size_t ifv = 0; ifv < m_nVectorFields; ++ifv) {
2830 recon.vector->at(
id, ifv) += modes[ifv] * reconstructionCoeffs[m_nScalarFields+ifv][ir];
2835 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2836 m_modes[ir].clear();
2841 sum(&recon, m_mean);
2858void POD::buildFields(
const std::vector<std::vector<double>> &reconstructionCoeffs, PiercedStorage<double> &fields,
2859 const std::vector<std::size_t> &scalarIds,
const std::vector<std::size_t> &podscalarIds,
2860 const std::vector<std::array<std::size_t, 3>> &vectorIds,
const std::vector<std::size_t> &podvectorIds,
2861 const std::unordered_set<long> *targetCells)
2863 std::size_t nsf = scalarIds.size();
2864 std::size_t nvf = vectorIds.size();
2865 if (nsf == 0 && nvf == 0) {
2869 std::unordered_set<long> targetCellsStorage;
2871 for (
const Cell &cell : m_podkernel->getMesh()->getCells()) {
2872 targetCellsStorage.insert(cell.getId());
2874 targetCells = &targetCellsStorage;
2878 for (
const long id : *targetCells) {
2879 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2880 double *recon = fields.rawData(rawIndex);
2881 for (std::size_t ifs = 0; ifs < m_nScalarFields; ++ifs) {
2882 double *reconsi = recon + scalarIds[ifs];
2885 for (std::size_t ifv = 0; ifv < m_nVectorFields; ++ifv) {
2886 for (std::size_t j = 0; j < 3; ++j) {
2887 double *reconvi = recon + vectorIds[ifv][j];
2894 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
2895 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2899 for (
const long id : *targetCells) {
2900 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
2901 double *recon = fields.rawData(rawIndex);
2903 if (m_nScalarFields > 0) {
2904 const double *modes = m_modes[ir].scalar->rawData(rawIndex);
2905 for (std::size_t ifs = 0; ifs < m_nScalarFields; ++ifs) {
2906 double *reconsi = recon + scalarIds[ifs];
2907 const double *modesi = modes + podscalarIds[ifs];
2908 (*reconsi) += (*modesi) * reconstructionCoeffs[ifs][ir];
2912 if (m_nVectorFields > 0) {
2913 const std::array<double,3> *modev = m_modes[ir].vector->rawData(rawIndex);
2914 for (std::size_t ifv = 0; ifv < m_nVectorFields; ++ifv) {
2915 const std::array<double,3> * modevi = modev + podvectorIds[ifv];
2916 for (std::size_t j = 0; j < 3; ++j) {
2917 double *reconvi = recon + vectorIds[ifv][j];
2918 (*reconvi) += (*modevi)[j] * reconstructionCoeffs[m_nScalarFields + ifv][ir];
2924 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2925 m_modes[ir].clear();
2930 sum(fields, m_mean, scalarIds, podscalarIds, vectorIds, podvectorIds, targetCells);
2941 throw std::runtime_error(
"POD: compute mapper can be called only in expert mode");
2942 _computeMapper(mesh);
2954 throw std::runtime_error(
"POD: update mapper can be called only in expert mode");
2955 _adaptionAlter(info);
2964 m_podkernel->computeMapper(mesh);
2965 m_podkernel->setMapperDirty(!m_expert);
2975 m_podkernel->adaptionPrepare(info);
2983void POD::_adaptionAlter(
const std::vector<adaption::Info> & info)
2985 m_podkernel->adaptionAlter(info);
2994 m_podkernel->adaptionCleanUp(info);
3004 const std::vector<std::size_t> &scalarIds,
const std::vector<std::size_t> &podscalarIds,
3005 const std::vector<std::array<std::size_t, 3>> &vectorIds,
const std::vector<std::size_t> &podvectorIds,
3006 const std::unordered_set<long> *targetCells)
3010 std::size_t nsf = scalarIds.size();
3011 std::size_t nvf = vectorIds.size();
3013 for (
long id : *targetCells) {
3014 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
3015 double *datag = fields.
rawData(rawIndex);
3016 double *datams = mode.
scalar->rawData(rawIndex);
3017 for (std::size_t i = 0; i < nsf; ++i) {
3018 double *datagi = datag + scalarIds[i];
3019 double *datamsi = datams + podscalarIds[i];
3020 (*datagi) -= (*datamsi);
3022 std::array<double,3>* datamv = mode.
vector->rawData(rawIndex);
3023 for (std::size_t i = 0; i < nvf; ++i) {
3024 std::array<double,3>* datamvi = datamv + podvectorIds[i];
3025 for (std::size_t j = 0; j < 3; ++j) {
3026 double *datagi = datag + vectorIds[i][j];
3027 (*datagi) -= (*datamvi)[j];
3041void POD::sum(PiercedStorage<double> &fields,
const pod::PODMode &mode,
3042 const std::vector<std::size_t> &scalarIds,
const std::vector<std::size_t> &podscalarIds,
3043 const std::vector<std::array<std::size_t, 3>> &vectorIds,
const std::vector<std::size_t> &podvectorIds,
3044 const std::unordered_set<long> *targetCells)
3048 std::size_t nsf = scalarIds.size();
3049 std::size_t nvf = vectorIds.size();
3051 for (
long id : *targetCells) {
3052 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(
id);
3053 double *datag = fields.rawData(rawIndex);
3054 double *datams = mode.scalar->rawData(rawIndex);
3055 for (std::size_t i = 0; i < nsf; ++i) {
3056 double *datagi = datag + scalarIds[i];
3057 double *datamsi = datams + podscalarIds[i];
3058 (*datagi) += (*datamsi);
3060 std::array<double,3>* datamv = mode.vector->rawData(rawIndex);
3061 for (std::size_t i = 0; i < nvf; ++i) {
3062 std::array<double,3>* datamvi = datamv + podvectorIds[i];
3063 for (std::size_t j = 0; j < 3; ++j) {
3064 double *datagi = datag + vectorIds[i][j];
3065 (*datagi) += (*datamvi)[j];
3081 log::cout() <<
"Writing snapshot ... " << std::endl;
3084 m_podkernel->getMesh()->getVTK().setDirectory(m_directory);
3087 std::vector<std::string > datastring;
3091 std::vector<std::vector<double>> fields(m_nScalarFields, std::vector<double>(field.
mesh->
getInternalCellCount(), 0.0));
3092 std::vector<std::vector<std::array<double,3>>> fieldv(m_nVectorFields, std::vector<std::array<double,3>>(field.
mesh->
getInternalCellCount(), {{0.0, 0.0, 0.0}}));
3093 std::vector<uint8_t> field_mask(field.mesh->getInternalCellCount(), 0);
3096 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
3098 for (
const Cell &cell : field.mesh->getCells()) {
3099 if (cell.isInterior()){
3100 long id = cell.getId();
3101 fields[isf][i] = field.scalar->at(
id, isf);
3106 field.mesh->getVTK().addData(m_nameScalarFields[isf], VTKFieldType::SCALAR, VTKLocation::CELL, fields[isf]);
3107 datastring.push_back(m_nameScalarFields[isf]);
3111 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
3113 for (
const Cell &cell : field.mesh->getCells()) {
3114 if (cell.isInterior()) {
3115 long id = cell.getId();
3116 fieldv[ivf][i] = field.vector->at(
id, ivf);
3120 std::string vname= m_nameVectorFields[ivf][0].substr(0,m_nameVectorFields[ivf][0].size()-2);
3122 field.mesh->getVTK().addData(vname, VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[ivf]);
3123 datastring.push_back(vname);
3128 for (
const Cell &cell : field.mesh->getCells()) {
3129 if (cell.isInterior()) {
3130 long id = cell.getId();
3131 field_mask[j] = m_filter[id];
3135 field.mesh->getVTK().addData(
"filter", VTKFieldType::SCALAR, VTKLocation::CELL, field_mask);
3136 datastring.push_back(
"filter");
3139 field.mesh->write();
3142 for (std::size_t is=0; is < datastring.size(); ++is) {
3143 field.mesh->getVTK().removeData(datastring[is]);
3155 log::cout() <<
"Writing mode ... " << std::endl;
3159 mode_mesh = m_podkernel->getMesh();
3163 std::vector<std::string> datastring;
3166 std::vector<std::vector<double>> fields(m_nScalarFields, std::vector<double>(mode_mesh->
getInternalCellCount(), 0.0));
3167 std::vector<std::vector<std::array<double,3>>> fieldv(m_nVectorFields, std::vector<std::array<double,3>>(mode_mesh->
getInternalCellCount(), {{0.0, 0.0, 0.0}}));
3168 std::vector<uint8_t> mask(mode_mesh->getInternalCellCount(), 0);
3171 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
3173 for (
const Cell &cell : mode_mesh->getCells()) {
3174 if (cell.isInterior()){
3175 long id = cell.getId();
3176 fields[isf][i] = m_modes[mode_index].scalar->at(
id, isf);
3181 mode_mesh->getVTK().addData(m_nameScalarFields[isf], VTKFieldType::SCALAR, VTKLocation::CELL, fields[isf]);
3182 datastring.push_back(m_nameScalarFields[isf]);
3186 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
3188 for (
const Cell &cell : mode_mesh->getCells()) {
3189 if (cell.isInterior()) {
3190 long id = cell.getId();
3191 fieldv[ivf][i] = m_modes[mode_index].vector->at(
id, ivf);
3195 std::string vname= m_nameVectorFields[ivf][0].substr(0,m_nameVectorFields[ivf][0].size()-2);
3197 mode_mesh->getVTK().addData(vname, VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[ivf]);
3198 datastring.push_back(vname);
3203 for (
const Cell &cell : mode_mesh->getCells()) {
3204 if (cell.isInterior()) {
3205 long id = cell.getId();
3206 mask[j] = m_filter[id];
3210 mode_mesh->getVTK().addData(
"mask", VTKFieldType::SCALAR, VTKLocation::CELL, mask);
3211 datastring.push_back(
"mask");
3217 for (std::size_t is=0; is < datastring.size(); ++is) {
3218 mode_mesh->getVTK().removeData(datastring[is]);
The Cell class defines the cells.
std::istream & getStream()
std::ostream & getStream()
The PODVolOctree is the specialized class of PODKernel for VolOctree meshes.
void setErrorThreshold(double threshold)
void computeMapper(const VolumeKernel *mesh)
std::size_t getModeCount()
void adaptionAlter(const std::vector< adaption::Info > &info)
void fillListActiveIDs(const PiercedStorage< bool > &bfield)
std::vector< std::array< std::string, 3 > > getVectorNames()
std::size_t getSnapshotCount()
void setTargetErrorFields(const std::vector< std::string > &namesf, const std::vector< std::array< std::string, 3 > > &namevf)
void evalReconstruction()
void unsetLeave1outSnapshots()
const std::string & getDirectory()
void dumpField(const std::string &name, const pod::PODField &field) const
const std::vector< pod::PODMode > & getModes()
void setRunMode(RunMode mode)
const VolumeKernel * getMesh()
std::vector< double > fieldsMax(pod::PODField &snap)
std::unique_ptr< PODKernel > & getKernel()
std::vector< std::string > getScalarNames()
const std::unordered_set< long int > & getListActiveIDs()
void setMeshType(MeshType type)
void setStaticMesh(bool flag)
void setSensorMask(const PiercedStorage< bool > &mask, VolumeKernel *mesh=nullptr)
WriteMode
Output Write Mode of the POD object. It defines the amount of information written by the POD object.
std::vector< double > fieldsl2norm(pod::PODField &snap)
void setMemoryMode(MemoryMode mode)
void write(const pod::PODField &snap, std::string file_name) const
double getErrorThreshold()
void setDirectory(const std::string &directory)
void adaptionPrepare(const std::vector< adaption::Info > &info)
void adaptionCleanUp(const std::vector< adaption::Info > &info)
const std::string & getName()
POD(MPI_Comm comm=MPI_COMM_WORLD)
void restoreDecomposition()
void readSnapshot(const pod::SnapshotFile &snap, pod::PODField &fieldr)
ReconstructionMode getReconstructionMode()
std::vector< std::string > getFieldsNames()
void reconstructFields(pod::PODField &field, pod::PODField &recon)
void evalErrorBoundingBox()
MemoryMode getMemoryMode()
void setExpert(bool mode=true)
void setUseMean(bool flag)
void addReconstructionSnapshot(const std::string &directory, const std::string &name)
void removeLeave1outSnapshot(const std::string &directory, const std::string &name)
std::vector< std::vector< double > > projectField(pod::PODField &field)
MeshType
Type of the Mesh used to compute the POD basis.
ReconstructionMode
Mode of Reconstruction of fields by using the POD basis.
const pod::PODMode & getMean()
MemoryMode
Memory Mode of the POD object. It defines the use of the memory resources.
void setSnapshots(const std::vector< pod::SnapshotFile > &database)
void setWriteMode(WriteMode mode)
std::size_t getListIDInternalCount()
void setReconstructionMode(ReconstructionMode mode)
ErrorMode
Mode of Error evaluation of a reconstructed fields by the POD basis.
void setEnergyLevel(double energy)
void setModeCount(std::size_t nmodes)
void setName(const std::string &name)
void setErrorMode(ErrorMode mode)
void setMesh(const std::string &directory, const std::string &name)
std::vector< std::vector< double > > getReconstructionCoeffs()
void addSnapshot(const std::string &directory, const std::string &name)
RunMode
Run Mode of the POD object. It defines if the POD basis has to be computed or restored.
PiercedVector< Cell > & getCells()
VTKUnstructuredGrid & getVTK()
void setVTKWriteTarget(WriteTarget targetCells)
virtual long getCellCount() const
bool dump(std::ostream &stream)
void write(VTKWriteMode mode=VTKWriteMode::DEFAULT)
long getInternalCellCount() const
Metafunction for generating a pierced storage.
__PS_POINTER__ rawData(std::size_t pos, std::size_t offset=0)
void setDirectory(const std::string &)
VTKField & addData(VTKField &&field)
void removeData(const std::string &)
void setName(const std::string &)
The VolOctree defines a Octree patch.
The VolumeKernel class provides an interface for defining volume patches.
void sum(const std::array< T, d > &x, T1 &s)
double norm(const std::array< T, d > &x, int p)
T dotProduct(const std::array< T, d > &x, const std::array< T, d > &y)
std::array< T, d > max(const std::array< T, d > &x, const std::array< T, d > &y)
std::array< T, d > min(const std::array< T, d > &x, const std::array< T, d > &y)
void write(std::ostream &stream, const std::vector< bool > &container)
void read(std::istream &stream, std::vector< bool > &container)
Logger & cout(log::Level defaultSeverity, log::Visibility defaultVisibility)
Logger & info(log::Visibility defaultVisibility)
PiercedStorage< std::array< double, 3 > > VectorStorage
PiercedStorage< double > ScalarStorage
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
void setMesh(VolumeKernel *lmesh)
The PODMode structure is used to store the modes inside pod classes.
std::unique_ptr< pod::VectorStorage > vector
std::unique_ptr< pod::ScalarStorage > scalar
The SnapFile structure is used to store the file names inside POD classes.