Loading...
Searching...
No Matches
pod.cpp
1/*---------------------------------------------------------------------------*\
2 *
3 * bitpit
4 *
5 * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6 *
7 * -------------------------------------------------------------------------
8 * License
9 * This file is part of bitpit.
10 *
11 * bitpit is free software: you can redistribute it and/or modify it
12 * under the terms of the GNU Lesser General Public License v3 (LGPL)
13 * as published by the Free Software Foundation.
14 *
15 * bitpit is distributed in the hope that it will be useful, but WITHOUT
16 * ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
17 * FITNESS FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public
18 * License for more details.
19 *
20 * You should have received a copy of the GNU Lesser General Public License
21 * along with bitpit. If not, see <http://www.gnu.org/licenses/>.
22 *
23\*---------------------------------------------------------------------------*/
24
25#include <cassert>
26#include <sstream>
27#include <sys/types.h>
28#include <sys/stat.h>
29#include <typeinfo>
30#include <unordered_map>
31#include <unordered_set>
32
33#if BITPIT_ENABLE_MPI
34# include <mpi.h>
35#endif
36
37#include "bitpit_private_lapacke.hpp"
38#include "bitpit_voloctree.hpp"
39
40#include "pod.hpp"
41#include "pod_voloctree.hpp"
42
43namespace bitpit {
44
59
63# if BITPIT_ENABLE_MPI
67POD::POD(MPI_Comm comm)
68# else
70# endif
71: m_filter(1), m_sensorMask(1)
72{
73
74#if BITPIT_ENABLE_MPI
75 m_communicator = MPI_COMM_NULL;
76#endif
77
78 m_podkernel = nullptr;
79 m_meshType = MeshType::UNDEFINED;
80 m_staticMesh = true; //default static mesh TODO change it
81 m_useMean = true; //default use mean fields
82 m_toUpdate = false; //default, switched to true during set<...>
83 m_nSnapshots = 0;
84 m_nModes = std::numeric_limits<std::size_t>::max();
85 m_nFields = 0;
86 m_nScalarFields = 0;
87 m_nVectorFields = 0;
88 m_nReconstructionSnapshots = 0;
89 m_sizeInternal = 0;
90 m_errorThreshold = 0;
91 setEnergyLevel(100.);
92
93# if BITPIT_ENABLE_MPI
94 initializeCommunicator(comm);
95 MPI_Comm_size(m_communicator, &m_nProcs);
96 MPI_Comm_rank(m_communicator, &m_rank);
97#else
98 m_rank = 0;
99 m_nProcs = 1;
100#endif
101
102 m_memoryMode = MemoryMode::MEMORY_NORMAL;
103 m_runMode = RunMode::COMPUTE;
104 m_writeMode = WriteMode::DUMP;
105 m_reconstructionMode = ReconstructionMode::MINIMIZATION;
106 m_errorMode = ErrorMode::NONE;
107 m_expert = false;
108
109 m_directory = ".";
110 m_name = "pod";
111
112}
113
118{
119 clear();
120}
121
126{
127
128 m_meshType = MeshType::UNDEFINED;
129 m_staticMesh = true;
130 m_useMean = true;
131
132 m_correlationMatrices.clear();
133 m_modes.clear();
134 m_minimizationMatrices.clear();
135
136 m_database.clear();
137 m_reconstructionDatabase.clear();
138 m_leave1outOffDatabase.clear();
139 m_toUpdate = false;
140 m_nSnapshots = 0;
141 m_nScalarFields = 0;
142 m_nVectorFields = 0;
143 m_nModes = 0;
144 m_nReconstructionSnapshots = 0;
145
146 if (m_podkernel) {
147 m_podkernel->clearMapper();
148 }
149
150# if BITPIT_ENABLE_MPI
151 freeCommunicator();
152# endif
153}
154
160void POD::setName(const std::string &name)
161{
162 m_name = name;
163}
164
170const std::string & POD::getName()
171{
172 return m_name;
173}
174
180void POD::setDirectory(const std::string &directory)
181{
182 struct stat info;
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!");
187 }
188
189 m_directory = directory + "/";
190}
191
197const std::string & POD::getDirectory()
198{
199 return m_directory;
200}
201
208void POD::addSnapshot(const std::string &directory, const std::string &name)
209{
210 pod::SnapshotFile file(directory, name);
211 addSnapshot(file);
212}
213
220{
221 m_database.push_back(file);
222 m_nSnapshots++;
223}
224
230void POD::setSnapshots(const std::vector<pod::SnapshotFile> &database)
231{
232 m_database = database;
233 m_nSnapshots = database.size();
234 m_toUpdate = true;
235}
236
243void POD::removeLeave1outSnapshot(const std::string &directory, const std::string &name)
244{
245 pod::SnapshotFile file(directory, name);
247}
248
255{
256 m_leave1outOffDatabase.push_back(file);
257}
258
263{
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);
271
272 }
273 }
274}
275
282void POD::addReconstructionSnapshot(const std::string &directory, const std::string &name)
283{
284 pod::SnapshotFile file(directory, name);
286}
287
294{
295 m_reconstructionDatabase.push_back(file);
296 m_nReconstructionSnapshots++;
297}
298
304void POD::setModeCount(std::size_t nmodes)
305{
306 m_nModes = nmodes;
307 m_toUpdate = true;
308}
309
315std::size_t POD::getModeCount()
316{
317 return m_nModes;
318}
319
326void POD::setEnergyLevel(double energy)
327{
328 m_energyLevel = energy - ENERGY_CHECK_TOLERANCE;
329}
330
337{
338 return m_energyLevel + ENERGY_CHECK_TOLERANCE;
339}
340
346void POD::setErrorThreshold(double threshold)
347{
348 m_errorThreshold = threshold;
349}
350
357{
358 return m_errorThreshold;
359}
360
367void POD::setTargetErrorFields(const std::vector<std::string> &namesf, const std::vector<std::array<std::string,3>> &namevf)
368{
369 std::size_t nsf= namesf.size();
370 std::size_t nvf= namevf.size();
371
372 m_nameTargetErrorFields.clear();
373
374 for (std::size_t ifield=0; ifield<nsf; ++ifield)
375 m_nameTargetErrorFields[namesf[ifield]] = ifield;
376
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;
380 }
381}
382
390{
391 if (m_meshType != MeshType::UNDEFINED)
392 throw std::runtime_error ("POD mesh type already set. Change not allowed.");
393
394 m_meshType = type;
395
396 switch (type)
397 {
399#if BITPIT_ENABLE_MPI
400 m_podkernel = std::unique_ptr<PODKernel>(new PODVolOctree(m_communicator));
401#else
402 m_podkernel = std::unique_ptr<PODKernel>(new PODVolOctree());
403#endif
404 break;
405
406 default:
407 throw std::runtime_error ("POD mesh type not allowed.");
408 break;
409 }
410}
411
418void POD::setMesh(std::unique_ptr<VolumeKernel> &&mesh)
419{
420 if (m_meshType == MeshType::UNDEFINED)
421 throw std::runtime_error ("POD mesh type not set.");
422
423 const VolOctree* _octreecast_mesh = dynamic_cast<const VolOctree*>(mesh.get());
424 if (_octreecast_mesh){
425 if (m_meshType != MeshType::VOLOCTREE)
426 throw std::runtime_error ("POD mesh type not set to VolOctree.");
427 m_podkernel->setMesh(std::move(mesh));
428 }
429 else{
430 throw std::runtime_error ("POD mesh type not allowed.");
431 }
432}
433
442void POD::setMesh(const std::string &directory, const std::string &name)
443{
444 pod::SnapshotFile file(directory, name);
445 setMesh(file);
446}
447
456{
457 if (m_meshType == MeshType::UNDEFINED)
458 throw std::runtime_error ("POD mesh type not set.");
459
460 std::unique_ptr<VolumeKernel> mesh = m_podkernel->readMesh(file);
461 setMesh(std::move(mesh));
462}
463
470{
471 return m_meshType;
472}
473
482void POD::setStaticMesh(bool flag)
483{
484 m_staticMesh = flag;
485}
486
492void POD::setUseMean(bool flag)
493{
494 m_useMean = flag;
495}
496
505{
506 if (m_memoryMode == mode) {
507 return;
508 }
509
510 switch (mode) {
511
513
514 m_memoryMode = MemoryMode::MEMORY_NORMAL;
515 break;
516
518
519 m_memoryMode = MemoryMode::MEMORY_LIGHT;
520 m_modes.clear();
521 break;
522
523 }
524}
525
530{
531 return m_memoryMode;
532}
533
542{
543 if (m_runMode == mode) {
544 return;
545 }
546
547 switch (mode) {
548
549 case RunMode::COMPUTE:
550
551 m_runMode = RunMode::COMPUTE;
552 break;
553
554 case RunMode::RESTORE:
555
556 m_runMode = RunMode::RESTORE;
557 break;
558
559 }
560}
561
566{
567 return m_runMode;
568}
569
576{
577 if (m_writeMode == mode)
578 return;
579
580 m_writeMode = mode;
581}
582
587{
588 return m_writeMode;
589}
590
598{
599 if (m_reconstructionMode == mode)
600 return;
601
602 switch (mode) {
603
605
606 m_reconstructionMode = ReconstructionMode::MINIMIZATION;
607 break;
608
610
611 m_reconstructionMode = ReconstructionMode::PROJECTION;
612 break;
613
614 }
615}
616
623{
624 return m_reconstructionMode;
625}
626
635{
636 if (m_errorMode == mode)
637 return;
638
639 m_errorMode = mode;
640}
641
648{
649 return m_errorMode;
650}
651
659void POD::setExpert(bool mode)
660{
661 m_expert = mode;
662}
663
672{
673 m_sensorMask.unsetKernel();
674 m_sensorMask.setStaticKernel(&(m_podkernel->getMesh()->getCells()));
675
676 if (m_staticMesh || mesh == m_podkernel->getMesh()){
677
678 for (const Cell &cell : m_podkernel->getMesh()->getCells()) {
679 long id = cell.getId();
680 m_sensorMask[id] = mask[id];
681 }
682
683 }
684 else{
685
686 if (mesh == nullptr)
687 throw std::runtime_error ("POD: null mesh pointer passed to setSensorMask");
688
689 //If mapping dirty or not computed, compute mapping of input mesh on pod mesh
690 if (m_podkernel->isMapperDirty())
691 _computeMapper(mesh);
692
693 std::unordered_set<long> trueCells;
694 for (const Cell & cell : mesh->getCells()){
695 long id = cell.getId();
696 if (mask[id])
697 trueCells.insert(id);
698 }
699
700 //Map true cells
701 std::unordered_set<long> mappedCells = m_podkernel->mapCellsToPOD(&trueCells);
702
703 m_podkernel->mapBoolFieldToPOD(mask, mesh, &mappedCells, m_sensorMask);
704
705 }
706
707 fillListActiveIDs(m_sensorMask);
708 m_toUpdate = true;
709
710}
711
718{
719 return m_nSnapshots;
720}
721
727std::vector<std::string> POD::getScalarNames()
728{
729 return m_nameScalarFields;
730}
731
737std::vector<std::array<std::string,3>> POD::getVectorNames()
738{
739 return m_nameVectorFields;
740}
741
748std::vector<std::string> POD::getFieldsNames()
749{
750 std::vector<std::string> names;
751 for (const std::string &ss : m_nameScalarFields)
752 names.push_back(ss);
753 for (const std::array<std::string,3> &ssa : m_nameVectorFields){
754 for (const std::string &ss : ssa)
755 names.push_back(ss);
756 }
757
758 return names;
759}
760
767{
768 if (m_podkernel == nullptr)
769 throw std::runtime_error ("POD mesh not built");
770
771 return m_podkernel->getMesh();
772}
773
780{
781 return m_mean;
782}
783
789const std::vector<pod::PODMode> & POD::getModes()
790{
791 return m_modes;
792}
793
799std::vector<std::vector<double> > POD::getReconstructionCoeffs()
800{
801 return m_reconstructionCoeffs;
802}
803
807const std::unordered_set<long> & POD::getListActiveIDs()
808{
809 return m_listActiveIDs;
810}
811
816{
817 return m_sizeInternal;
818}
819
823std::unique_ptr<PODKernel> & POD::getKernel()
824{
825 return m_podkernel;
826}
827
836{
837 // Decomposition
838 if (m_runMode == RunMode::COMPUTE) {
839 log::cout() << "pod : evaluate decomposition ... " << std::endl;
841 } else if (m_runMode == RunMode::RESTORE) {
842 log::cout() << "pod : restore decomposition ... " << std::endl;
844 }
845
846 // Reconstruction
847 log::cout() << "pod : evaluate reconstruction ... " << std::endl;
849}
850
855{
856 // Dump decomposition
858}
859
864{
865 // Restore decomposition
867}
868
873{
874 // Evaluate mean field and mesh
875 log::cout() << "pod : computing mean field and pod mesh... " << std::endl;
876 evalMeanMesh();
877
878 // Initialize ID list
879 fillListActiveIDs(m_filter);
880
881 // Evaluate correlation matrices
882 log::cout() << "pod : computing correlation matrix... " << std::endl;
884
885 // Compute eigenvalues and eigenvectors
886 log::cout() << "pod : computing eigenvectors... " << std::endl;
887 evalEigen();
888
889 // Compute POD modes
890 log::cout() << "pod : computing pod modes... " << std::endl;
891 evalModes();
892
893 // Dump pod
894 if (m_writeMode != WriteMode::NONE) {
895 log::cout() << "pod : dumping... " << std::endl;
896 dump();
897 }
898}
899
904{
905 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
906
907 log::cout() << "Dumping POD ... " << std::endl;
908 // Dump the POD info (no reconstruction info)
909 std::string header = "pod";
910 OBinaryArchive binaryWriter(m_directory + m_name, ARCHIVE_VERSION, header, dumpBlock);
911 std::ostream &dataStream = binaryWriter.getStream();
912 utils::binary::write(dataStream, m_meshType);
913 utils::binary::write(dataStream, m_nModes);
914 utils::binary::write(dataStream, m_nSnapshots);
915 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
916 utils::binary::write(dataStream, m_database[i].directory);
917 utils::binary::write(dataStream, m_database[i].name);
918 }
919 utils::binary::write(dataStream,m_staticMesh);
920 utils::binary::write(dataStream,m_memoryMode);
921 binaryWriter.close();
922
923 // Dump the modes, mean and meshPOD
924 if (m_memoryMode == MemoryMode::MEMORY_NORMAL) {
925 for (std::size_t ir = 0; ir < m_nModes; ++ir)
926 dumpMode(ir);
927
928 }
929
930 // Dumping mean as snapshot
931 if (m_useMean){
932 dumpMode(m_nModes);
933 }
934
935 // Dumping mask
936 {
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();
942 }
943
944 // Dumping mesh
945 {
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();
950 }
951
952 // Write debug info
953 if (m_writeMode == WriteMode::DEBUG) {
954 std::vector<std::string > datastring;
955 m_podkernel->getMesh()->getVTK().setDirectory(m_directory);
956 m_podkernel->getMesh()->getVTK().setName(m_name);
957
958 std::size_t nf=m_nModes;
959 if (m_useMean)
960 nf=nf+1;
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}})));
963
964 if (m_useMean){
965 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
966 int i = 0;
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);
971 ++i;
972 }
973 }
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");
976 }
977 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
978 int i = 0;
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);
983 ++i;
984 }
985 }
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");
988 }
989 }
990 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
991 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
992 readMode(ir);
993
994 std::size_t nir = ir;
995 if (m_useMean)
996 nir=ir+1;
997
998 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
999 int i=0;
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);
1004 ++i;
1005 }
1006 }
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));
1009 }
1010
1011 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
1012 int i=0;
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);
1017 ++i;
1018 }
1019 }
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));
1022 }
1023
1024 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
1025 m_modes[ir].clear();
1026 }
1027
1028 std::vector<uint8_t> mask(m_podkernel->getMesh()->getInternalCellCount(), 0);
1029 int i=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];
1034 ++i;
1035 }
1036 }
1037
1038 m_podkernel->getMesh()->getVTK().addData("filter", VTKFieldType::SCALAR, VTKLocation::CELL, mask);
1039 datastring.push_back("filter");
1040 m_podkernel->getMesh()->write();
1041
1042 for (std::size_t is = 0; is < datastring.size(); ++is)
1043 m_podkernel->getMesh()->getVTK().removeData(datastring[is]);
1044 }
1045}
1046
1051{
1052 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
1053
1054 // Restore the POD info (no reconstruction info)
1055 log::cout() << "Restore POD ... " << std::endl;
1056
1057 std::string header = "pod";
1058 IBinaryArchive binaryReader(m_directory + m_name, dumpBlock);
1059 std::istream &dataStream = binaryReader.getStream();
1060 MeshType m_meshType_;
1061 utils::binary::read(dataStream, m_meshType_);
1062 setMeshType(m_meshType_);
1063 std::size_t nr_;
1064 utils::binary::read(dataStream, nr_);
1065 m_nModes = std::min(m_nModes, nr_);
1066 m_modes.resize(m_nModes);
1067 utils::binary::read(dataStream, m_nSnapshots);
1068 m_database.resize(m_nSnapshots);
1069 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
1070 utils::binary::read(dataStream, m_database[i].directory);
1071 utils::binary::read(dataStream, m_database[i].name);
1072 }
1073 utils::binary::read(dataStream,m_staticMesh);
1074 utils::binary::read(dataStream,m_memoryMode);
1075 binaryReader.close();
1076
1077 // Restore mesh
1078 {
1079 pod::SnapshotFile podfile;
1080 podfile.directory = m_directory;
1081 podfile.name = m_name;
1082 m_podkernel->restoreMesh(podfile);
1083 }
1084
1085 // Restore the modes, mean and mesh
1086 {
1087 m_filter.unsetKernel();
1088 m_filter.setStaticKernel(&m_podkernel->getMesh()->getCells());
1089 m_sensorMask.unsetKernel();
1090 m_sensorMask.setStaticKernel(&m_podkernel->getMesh()->getCells());
1091
1092 if (m_memoryMode == MemoryMode::MEMORY_NORMAL) {
1093 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
1094 readMode(ir);
1095 }
1096 }
1097 }
1098
1099 // Restore mean
1100 if (m_useMean){
1101 readMode(m_nModes);
1102 }
1103
1104 // Restore mask
1105 {
1106 IBinaryArchive binaryReader(m_directory + m_name + ".mask.data", dumpBlock);
1107 std::istream &dataStream = binaryReader.getStream();
1108 m_sensorMask.restore(dataStream);
1109 binaryReader.close();
1110 }
1111
1112 // Lestore ID list
1113 fillListActiveIDs(m_sensorMask);
1114 m_toUpdate = true;
1115}
1116
1121{
1125
1126 // Evaluate mean field and mesh
1127 log::cout() << "pod : computing pod mesh... " << std::endl;
1128 evalMeanMesh();
1129
1130 // Initialize ID list
1131 fillListActiveIDs(m_filter);
1132
1133 // Evaluate correlation matrices
1134 log::cout() << "pod : computing correlation matrix... " << std::endl;
1135 evalCorrelation();
1136
1137 // Initialize error field
1138 initErrorMaps();
1139
1140 // Initialize help variables
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();
1145
1146 m_nSnapshots=h_nSnapshots-1;
1147
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;
1151
1152 m_reconstructionDatabase.clear();
1153 m_nReconstructionSnapshots = 0;
1154 addReconstructionSnapshot(m_database[id]);
1155 m_database.erase(m_database.begin()+id);
1156
1157 for (std::size_t ifield = 0; ifield<m_nFields; ++ifield){
1158 // get rid of id row
1159 m_correlationMatrices[ifield].erase(m_correlationMatrices[ifield].begin()+id*h_nSnapshots, m_correlationMatrices[ifield].begin()+id*h_nSnapshots+h_nSnapshots);
1160
1161 // get rid of id column
1162 std::size_t it=id;
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;
1166 }
1167 }
1168
1169 // Compute eigenvalues and eigenvectors
1170 log::cout() << "pod : computing eigenvectors... " << std::endl;
1171 evalEigen();
1172
1173 // Compute POD modes
1174 log::cout() << "pod : computing pod modes... " << std::endl;
1175 evalModes();
1176
1177 // Evaluate reconstruction & compute error fields
1179
1180 // Reassignment for future usage
1181 m_correlationMatrices=h_correlationMatrices;
1182 m_database=h_database;
1183 m_toUpdate=true;
1184 }
1185
1186 //dump error
1187 if (m_writeMode != WriteMode::NONE){
1188 std::string ename = m_name + ".error";
1189 dumpField(ename, m_errorMap);
1190 }
1191 m_nSnapshots=h_nSnapshots;
1192}
1193
1200{
1201 if (m_database.empty())
1202 throw std::runtime_error ("POD database empty");
1203
1204 // If the mesh is static do nothing on mesh and compute mean fields directly
1205 _evalMeanMesh();
1206
1207 // Set mesh POD to write only internal cells
1208#if BITPIT_ENABLE_MPI
1209 m_podkernel->getMesh()->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_INTERNAL);
1210#else
1211 m_podkernel->getMesh()->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_ALL);
1212#endif
1213}
1214
1218void POD::_evalMeanMesh()
1219{
1220
1221 if (m_podkernel->getMesh() == nullptr){
1222
1223 //Read first mesh and use as meshPOD
1224 std::unique_ptr<VolumeKernel> readMesh = m_podkernel->readMesh(m_database[0]);
1225 m_podkernel->setMesh(std::move(readMesh));
1226
1227 //Dynamic mesh case
1228 if (!m_staticMesh){
1229 // Compute meshPOD (starting from inital mesh) and fill filter
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();
1235 }
1236 }
1237 }
1238
1239 {
1240 // Read first snapshot to set n scalar and vector fields
1241 pod::PODField readf;
1242 readSnapshot(m_database[0], readf);
1243 }
1244
1245 //Compute cells volume
1246 m_podkernel->evalCellsVolume();
1247
1248 m_filter.unsetKernel();
1249 m_filter.setStaticKernel(&(m_podkernel->getMesh()->getCells()));
1250 m_filter.fill(true);
1251
1252
1253 if (m_useMean){
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}});
1258
1259 // Compute mean (starting from first snapshot) & filter+mask
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);
1264 if (!m_staticMesh){
1265 //Compute mapping of read mesh on pod mesh
1266 _computeMapper(readf.mesh);
1267 readf = pod::PODField(m_podkernel->mapPODFieldToPOD(readf, nullptr));
1268 m_podkernel->clearMapper();
1269 }
1270
1271 for (const Cell &cell : m_podkernel->getMesh()->getCells()) {
1272 m_filter[cell.getId()] = (m_filter[cell.getId()]) && (readf.mask->at(cell.getId())); //N.B. mask == filter
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);
1278 }
1279 else{
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}};
1284 }
1285 }
1286 }
1287 }
1288 else {
1289 // filter+mask
1290 for (std::size_t i = 0; i < m_nSnapshots; ++i) {
1291 pod::PODField readf;
1292 readSnapshot(m_database[i], readf);
1293 if (!m_staticMesh)
1294 {
1295 _computeMapper(readf.mesh);
1296 readf = pod::PODField(m_podkernel->mapPODFieldToPOD(readf, nullptr));
1297 m_podkernel->clearMapper();
1298 }
1299
1300 for (const Cell &cell : m_podkernel->getMesh()->getCells())
1301 m_filter[cell.getId()] = (m_filter[cell.getId()]) && (readf.mask->at(cell.getId())); //N.B. mask == filter
1302
1303 }
1304 }
1305
1306 setSensorMask(m_filter, m_podkernel->getMesh());
1307}
1308
1316{
1317
1318 //NOTE!! NO GHOSTS IN ACTIVEIDS (*)
1319 m_listActiveIDs.clear();
1320 std::unordered_set<long> listGhost;
1321 for (const Cell &cell : m_podkernel->getMesh()->getCells()) {
1322 long id = cell.getId();
1323 //TODO reset to commented when solid/fluid field on ghosts is dumped
1324 //if (mask[id]) {
1325 if (mask[id] && cell.isInterior()) {
1326 if (cell.isInterior()) {
1327 m_listActiveIDs.insert(id);
1328 }else{
1329 listGhost.insert(id);
1330 }
1331 }
1332 }
1333
1334 m_sizeInternal = m_listActiveIDs.size();
1335
1336 // (*)
1337 //m_listActiveIDs.insert(listGhost.begin(), listGhost.end());
1338
1339}
1340
1345{
1346 initCorrelation();
1347 for (std::size_t i = 0; i < m_nSnapshots; ++i){
1348 pod::PODField snapi;
1349 readSnapshot(m_database[i], snapi);
1350 if (!m_staticMesh){
1351 _computeMapper(snapi.mesh);
1352 snapi = pod::PODField(m_podkernel->mapPODFieldToPOD(snapi, nullptr));
1353 m_podkernel->clearMapper();
1354 }
1355
1356 if (m_useMean)
1357 diff(&snapi, m_mean);
1358
1359 for (std::size_t j = i; j < m_nSnapshots; ++j){
1360 log::cout() << "pod : evaluation " << i << "," << j << " term of correlation matrix " << std::endl;
1361 pod::PODField snapj;
1362 readSnapshot(m_database[j], snapj);
1363 if (!m_staticMesh){
1364 _computeMapper(snapj.mesh);
1365 snapj = pod::PODField(m_podkernel->mapPODFieldToPOD(snapj, nullptr));
1366 m_podkernel->clearMapper();
1367 }
1368
1369 if (m_useMean)
1370 diff(&snapj, m_mean);
1371
1372 evalCorrelationTerm(i, snapi, j, snapj);
1373
1374 }
1375 }
1376
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);
1380# endif
1381
1382}
1383
1387void POD::initCorrelation()
1388{
1389 m_correlationMatrices.clear();
1390 m_correlationMatrices.resize(m_nFields,std::vector<double>(m_nSnapshots*m_nSnapshots, 0.0));
1391}
1392
1396void POD::initErrorMaps()
1397{
1398 m_errorMap.clear();
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;
1408
1409}
1410
1414void POD::evalCorrelationTerm(int i, pod::PODField & snapi, int j, pod::PODField & snapj)
1415{
1416
1417 std::unordered_set<long>::iterator it;
1418 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it){
1419 long id = *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);
1426 datasi++;
1427 datasj++;
1428 }
1429 }
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);
1435 datavi++;
1436 datavj++;
1437 }
1438 }
1439
1440 }
1441
1442}
1443
1449{
1450 m_reconstructionCoeffs.clear();
1451 m_reconstructionCoeffs.resize(m_nFields, std::vector<double>(m_nModes,0.0));
1452
1453 for (std::size_t i = 0; i < m_nReconstructionSnapshots; ++i){
1454 pod::PODField snapi, reconi, erri;
1455 readSnapshot(m_reconstructionDatabase[i], snapi);
1456
1457 if (!m_staticMesh){
1458
1459 //Compute mapping of snapshot mesh on pod mesh
1460 _computeMapper(snapi.mesh);
1461
1462 // Map snapshot on pod mesh
1463 snapi = pod::PODField(m_podkernel->mapPODFieldToPOD(snapi, nullptr));
1464 m_podkernel->clearMapper();
1465
1466 }
1467
1468 reconstructFields(snapi, reconi);
1469 std::string rname = m_name + ".recon." + m_reconstructionDatabase[i].name;
1470
1471 if (m_writeMode != WriteMode::NONE && m_errorMode != ErrorMode::COMBINED)
1472 dumpField(rname, reconi);
1473
1474 if (m_errorMode != ErrorMode::COMBINED){
1475 // one coeff file for each reconstruction, each line is related to a field, each term is related to a mode
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;
1479
1480 outr.close();
1481 }
1482 if (m_errorMode == ErrorMode::SINGLE){
1483 std::vector<double> enorm = fieldsl2norm(m_errorMap);
1484 std::string ename = m_name + ".error."+ m_reconstructionDatabase[i].name;
1485
1486 if (m_writeMode != WriteMode::NONE)
1487 dumpField(ename, m_errorMap);
1488
1489 // one coeff file for each reconstruction, each line is related to a field
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;
1493
1494 oute.close();
1495 }
1496 }
1497}
1498
1503{
1505 pod::SnapshotFile efile("./pod", m_name+".error");
1506 log::cout()<< "Reading error field ..."<< std::endl;
1507 log::cout()<< efile.directory+"/"+efile.name +"..."<< std::endl;
1508 pod::PODField error;
1509 readSnapshot(efile, error);
1510
1511 fillListActiveIDs(*error.mask);
1512
1513 std::vector<std::size_t> scalarIds, vectorIds;
1514 std::map<std::string, std::size_t> targetErrorFields = m_nameTargetErrorFields;
1515
1516 // Find scalar target fields
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);
1523 }
1524 count++;
1525 }
1526
1527 // Find vector target fields
1528 for (const std::array<std::string,3> &valv : m_nameVectorFields) {
1529 std::size_t ic = 0;
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()) {
1534 toerase[ic] = val;
1535 ic++;
1536 }
1537
1538 if (ic == 3) {
1539 vectorIds.push_back(count);
1540 for (std::size_t i = 0; i < 3; ++i)
1541 targetErrorFields.erase(toerase[i]);
1542 }
1543 }
1544 count++;
1545 }
1546
1547 if (targetErrorFields.size() != 0) {
1548 std::string verror = " ";
1549 for (const auto &targetErrorField : targetErrorFields)
1550 verror += targetErrorField.first + " ";
1551
1552 throw std::runtime_error ("POD: fields not found in POD base or incomplete vector: " + verror);
1553 }
1554
1555 std::size_t nsf = scalarIds.size();
1556 std::size_t nvf = vectorIds.size();
1557
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}});
1561
1562 std::unordered_set<long>::iterator it;
1563 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it) {
1564 long id = *it;
1565 std::array<double, 3> cellCentroid = error.mesh->evalCellCentroid(id);
1566 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(id);
1567
1568 if (nsf) {
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]]);
1575 }
1576 }
1577 }
1578 if (nvf) {
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]]);
1585 }
1586 }
1587 }
1588 }
1589
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);
1593# endif
1594
1595 log::cout()<< ">> Fields Bounding boxes " << std::endl;
1596 log::cout()<< "min: "<< minBoxes << std::endl;
1597 log::cout()<< "max: "<< maxBoxes << std::endl;
1598
1599 std::array<double, 3> minBox = {{0.0, 0.0, 0.0}};
1600 std::array<double, 3> maxBox = {{0.0, 0.0, 0.0}};
1601
1602 for (std::size_t i =0; i < m_nFields; ++i){
1603 maxBox=max(maxBoxes[i],maxBox);
1604 minBox=min(minBoxes[i],minBox);
1605 }
1606
1607 log::cout()<< ">> Bounding box " << std::endl;
1608 log::cout()<< "("<< minBox << ") ("<< maxBox << ")"<< std::endl;
1609
1610 bool runSolver;
1611#if BITPIT_ENABLE_MPI
1612 runSolver = (m_rank == 0);
1613#else
1614 runSolver = true;
1615#endif
1616 if (runSolver){
1617 std::ofstream outBB(m_directory + m_name + ".box.dat", std::ofstream::out);
1618 outBB << minBox << std::endl;
1619 outBB << maxBox << std::endl;
1620 outBB.close();
1621 }
1622
1623}
1624
1631{
1632
1633 //set the mesh of reconstructed field
1634 reconi.mesh = field.mesh;
1635 reconi.mask = std::unique_ptr<PiercedStorage<bool>>(new PiercedStorage<bool>(*field.mask));
1636 evalReconstructionCoeffs(field);
1637 buildFields(m_reconstructionCoeffs, reconi);
1638 if (m_errorMode != ErrorMode::NONE){
1639 if (m_errorMode == ErrorMode::SINGLE)
1640 initErrorMaps();
1641 buildErrorMaps(field,reconi);
1642 }
1643
1644}
1645
1651void POD::reconstructFields(const std::vector<std::vector<double>> &reconstructionCoeffs,
1652 pod::PODField & reconi)
1653{
1654 reconi.setMesh(m_podkernel->getMesh());
1655 buildFields(reconstructionCoeffs, reconi);
1656}
1657
1667 std::map<std::string, std::size_t> targetFields,
1668 const std::unordered_set<long> * targetCells)
1669{
1670 std::vector<std::size_t> podscalarIds;
1671 std::vector<std::size_t> podvectorIds;
1672
1673 std::vector<std::size_t> scalarIds;
1674 std::vector<std::array<std::size_t, 3>> vectorIds;
1675
1676 // Find scalar target fields
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);
1685 count++;
1686 }
1687 }
1688
1689 // Find vector target fields
1690 count = 0;
1691 for (const std::array<std::string,3> &valv : m_nameVectorFields) {
1692 std::size_t ic = 0;
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;
1700 toerase[ic] = val;
1701 ic++;
1702 }
1703
1704 if (ic == 3) {
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]);
1709 }
1710 }
1711 count++;
1712 }
1713
1714 if (targetFields.size() != 0) {
1715 std::string verror = " ";
1716 for (const auto &targetField : targetFields)
1717 verror += targetField.first + " ";
1718
1719 throw std::runtime_error ("POD: fields not found in POD base or incomplete vector: " + verror);
1720 }
1721
1722
1723 if (m_staticMesh){
1724
1725 evalReconstructionCoeffs(fields, scalarIds, podscalarIds, vectorIds, podvectorIds);
1726
1727 buildFields(m_reconstructionCoeffs, fields, scalarIds, podscalarIds, vectorIds, podvectorIds, targetCells);
1728
1729 }
1730 else{
1731
1732 //If mapping not computed, compute mapping of input mesh on pod mesh
1733 if (m_podkernel->isMapperDirty())
1734 _computeMapper(mesh);
1735
1736 // Map data fields on pod mesh
1737 PiercedStorage<double> mappedFields = m_podkernel->mapFieldsToPOD(fields, mesh, &m_listActiveIDs, scalarIds, vectorIds);
1738
1739 evalReconstructionCoeffs(mappedFields, scalarIds, podscalarIds, vectorIds, podvectorIds);
1740
1741 //Map target cells
1742 std::unordered_set<long> mappedCells;
1743 std::unordered_set<long>* ptr_mappedCells = nullptr;
1744 if (targetCells){
1745 mappedCells = m_podkernel->mapCellsToPOD(targetCells);
1746 ptr_mappedCells = &mappedCells;
1747 }
1748
1749 buildFields(m_reconstructionCoeffs, mappedFields, scalarIds, podscalarIds, vectorIds, podvectorIds, ptr_mappedCells);
1750
1751 m_podkernel->mapFieldsFromPOD(fields, mesh, targetCells, mappedFields, scalarIds, vectorIds);
1752
1753 }
1754
1755}
1756
1763std::vector<std::vector<double>> POD::projectField(pod::PODField &field)
1764{
1765 std::vector<std::vector<double>> reconstructionCoeffs(m_nFields, std::vector<double>(m_nModes,0.0));
1766 if (m_useMean) {
1767 diff(&field, m_mean);
1768 }
1769
1770 // Outer cycle over the POD modes
1771 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
1772 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
1773 readMode(ir);
1774 }
1775
1776 // Outer cycle over active cells
1777 for (long id : m_listActiveIDs) {
1778 std::size_t rawIndex = field.mesh->getCells().getRawIndex(id);
1779
1780 // Inner cycle over scalar fields
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);
1786 }
1787 }
1788
1789 // Inner cycle over vector fields
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);
1795 }
1796 }
1797 }
1798 }
1799
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);
1803 }
1804#endif
1805
1806 if (m_useMean) {
1807 sum(&field, m_mean);
1808 }
1809
1810 return reconstructionCoeffs;
1811}
1812
1819void POD::buildErrorMaps(pod::PODField & snap, pod::PODField & recon)
1820{
1821 std::vector<double> norm = fieldsMax(snap);
1822
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]);
1831 errors++;
1832 recons++;
1833 snaps++;
1834 }
1835 }
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]);
1843
1844 errorv++;
1845 reconv++;
1846 snapv++;
1847 }
1848 }
1849 }
1850}
1851
1855void POD::initMinimization()
1856{
1857 m_minimizationMatrices.clear();
1858 m_minimizationMatrices.resize(m_nFields,std::vector<double>(m_nModes*m_nModes, 0.0));
1859}
1860
1864void POD::evalMinimizationMatrices()
1865{
1866 if (m_toUpdate){
1867 initMinimization();
1868
1869
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;
1874 }
1875
1876 }else{
1877
1878 for (std::size_t ir = 0; ir < m_nModes; ++ir){
1879 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
1880 readMode(ir);
1881
1882 for (std::size_t jr = ir; jr < m_nModes; ++jr){
1883 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
1884 readMode(jr);
1885
1886 std::unordered_set<long>::iterator it;
1887 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it){
1888 long id = *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);
1895 datasi++;
1896 datasj++;
1897 }
1898
1899 if (ir != jr){
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]; // sym: just in case...
1902 }
1903 }
1904 }
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);
1910 datavi++;
1911 datavj++;
1912 }
1913
1914 if (ir != jr){
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]; // sym: just in case...
1917 }
1918 }
1919 }
1920 }
1921 }
1922
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);
1926# endif
1927 }
1928 m_toUpdate = false;
1929 }
1930
1931}
1932
1938void POD::solveMinimization(std::vector<std::vector<double> > & rhs)
1939{
1940 bool runSolver;
1941#if BITPIT_ENABLE_MPI
1942 runSolver = (m_rank == 0);
1943#else
1944 runSolver = true;
1945#endif
1946
1947 if (runSolver) {
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);
1951
1952 for (std::size_t j = 0; j < m_nModes*m_nModes; ++j)
1953 A[j] = m_minimizationMatrices[i][j];
1954
1955 lapack_int info = LAPACKE_dgesv(LAPACK_COL_MAJOR, m_nModes, 1, A.data(), m_nModes, ipiv.data(), rhs[i].data(), m_nModes); //lapacke
1956 if (info != 0) {
1957 log::cout() << "WARNING! algorithm convergence info " << info << std::endl;
1958 throw std::runtime_error("Unable to solve minimization problem.");
1959 }
1960 }
1961 }
1962 else{
1963 for (std::size_t i = 0; i < m_nFields; ++i) {
1964 rhs[i] = std::vector<double>(m_nModes, 0.0);
1965 }
1966 }
1967}
1968
1973{
1974 m_lambda.clear();
1975 m_lambda.resize(m_nFields,std::vector<double>(m_nSnapshots, 0.0));
1976
1977 m_podCoeffs.clear();
1978 m_podCoeffs.resize(m_nFields, std::vector<std::vector<double>>(m_nSnapshots, std::vector<double>(m_nSnapshots,0.0)));
1979
1980 bool runSolver;
1981#if BITPIT_ENABLE_MPI
1982 runSolver = (m_rank == 0);
1983#else
1984 runSolver = true;
1985#endif
1986
1987 if (runSolver) {
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];
1994
1995 lapack_int info = LAPACKE_dsyev(LAPACK_COL_MAJOR, 'V', 'U', m_nSnapshots, Marr.data(), m_nSnapshots, alambda.data()); //lapacke
1996 if (info != 0) {
1997 log::cout() << "WARNING! algorithm convergence info " << info << std::endl;
1998 throw std::runtime_error("Unable to solve minimization problem.");
1999 }
2000
2001 // Check number of modes with energy level
2002 checkModeCount(alambda.data(), i);
2003
2004 for (std::size_t n = 0; n < m_nSnapshots; ++n)
2005 m_lambda[i][n] = alambda[m_nSnapshots - n - 1];
2006
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]));
2010 }
2011 }
2012 if (m_errorMode != ErrorMode::COMBINED){
2013 // one coeff file for each mode, each line is related to a field, each term is related to a snapshot
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;
2018
2019 out.close();
2020 }
2021
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;
2025
2026 outl.close();
2027 }
2028 }
2029
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);
2037
2038 MPI_Bcast(m_lambda[i].data(), m_nSnapshots, MPI_DOUBLE, 0, m_communicator);
2039 }
2040#endif
2041}
2042
2049void POD::checkModeCount(double *alambda, std::size_t ifield)
2050{
2051 if (ifield == 0){
2052 _m_nr.clear();
2053 if (m_useMean)
2054 m_nModes = std::min(m_nModes, m_nSnapshots-1);
2055 else
2056 m_nModes = std::min(m_nModes, m_nSnapshots);
2057 }
2058
2059 std::size_t nr = 0;
2060 std::size_t hn = 0;
2061 if (!m_useMean)
2062 hn=1;
2063 double en = 0;
2064 double sumen = 0;
2065 std::vector<double> dlambda(m_nSnapshots);
2066
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]);
2070 }
2071 while (nr < m_nModes && en < m_energyLevel){
2072 en += std::abs(dlambda[nr]) / sumen * 100;
2073 nr++;
2074 }
2075 _m_nr.push_back(nr);
2076 if (ifield == m_nFields-1)
2077 m_nModes = *(std::max_element(_m_nr.begin(), _m_nr.end()));
2078}
2079
2084{
2085 m_modes.resize(m_nModes);
2086 _evalModes();
2087}
2088
2092void POD::_evalModes()
2093{
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}});
2099 }
2100
2101 for (std::size_t is = 0; is < m_nSnapshots; ++is){
2102 pod::PODField snapi;
2103 readSnapshot(m_database[is], snapi);
2104 if (!m_staticMesh){
2105 _computeMapper(snapi.mesh);
2106 snapi = pod::PODField(m_podkernel->mapPODFieldToPOD(snapi, nullptr));
2107 m_podkernel->clearMapper();
2108 }
2109
2110 if (m_useMean)
2111 diff(&snapi, m_mean);
2112
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];
2121 modes++;
2122 datas++;
2123 }
2124 }
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];
2130 modev++;
2131 datav++;
2132 }
2133 }
2134 }
2135 if (m_memoryMode == MemoryMode::MEMORY_LIGHT){
2136 dumpMode(ir);
2137 m_modes[ir].clear();
2138 }
2139 }
2140 }
2141}
2142
2150{
2151
2152 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2153
2154 std::string filename = std::string(snap.directory) + "/" + std::string(snap.name) + ".data";
2155 IBinaryArchive binaryReader(filename, dumpBlock);
2156 std::istream &dataStream = binaryReader.getStream();
2157
2158 //Set mesh field
2159 if (m_staticMesh){
2160 fieldr.setMesh(m_podkernel->getMesh());
2161 }
2162 else{
2163 //Read and set the mesh
2164 fieldr.setMesh(m_podkernel->readMesh(snap));
2165 }
2166 VolumeKernel * mesh = fieldr.mesh;
2167
2168 // Restore solved cells
2169 fieldr.mask = std::unique_ptr<PiercedStorage<bool>>(new PiercedStorage<bool>(1, &mesh->getCells()));
2170 fieldr.mask->restore(dataStream);
2171
2172 // For each snapshot the global number of fields are read & updated again .
2173
2174 // Restore scalar fields
2175 utils::binary::read(dataStream, m_nScalarFields);
2176 m_nameScalarFields.resize(m_nScalarFields);
2177 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2178 utils::binary::read(dataStream, m_nameScalarFields[i]);
2179
2180 if (m_nScalarFields) {
2181 fieldr.scalar = std::unique_ptr<pod::ScalarStorage>(new pod::ScalarStorage(m_nScalarFields, &mesh->getCells()));
2182 fieldr.scalar->restore(dataStream);
2183 }
2184
2185 // Restore vector fields
2186 utils::binary::read(dataStream, m_nVectorFields);
2187 m_nFields = m_nScalarFields + m_nVectorFields;
2188 m_nameVectorFields.resize(m_nVectorFields);
2189 for (std::size_t i = 0; i < m_nVectorFields; ++i)
2190 utils::binary::read(dataStream, m_nameVectorFields[i]);
2191
2192 if (m_nVectorFields) {
2193 fieldr.vector = std::unique_ptr<pod::VectorStorage>(new pod::VectorStorage(m_nVectorFields, &mesh->getCells()));
2194 fieldr.vector->restore(dataStream);
2195 }
2196 binaryReader.close();
2197
2198 // Write vtk with field in DEBUG mode
2199 if (m_writeMode == WriteMode::DEBUG) {
2200 std::vector<std::string > datastring;
2201 mesh->getVTK().setDirectory(m_directory);
2202 mesh->getVTK().setName(m_name + "." + snap.name);
2203 std::vector<uint8_t> mask(mesh->getCellCount());
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()));
2206
2207 int i=0;
2208 for (const Cell &cell : mesh->getCells()) {
2209 long id = cell.getId();
2210 mask[i] = uint8_t(fieldr.mask->at(id));
2211 ++i;
2212 }
2213 mesh->getVTK().addData("mask", VTKFieldType::SCALAR, VTKLocation::CELL, mask);
2214 datastring.push_back("mask");
2215
2216 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
2217 int i = 0;
2218 for (const Cell &cell : mesh->getCells()) {
2219 long id = cell.getId();
2220 fields[isf][i] = fieldr.scalar->at(id, isf);
2221 ++i;
2222 }
2223 mesh->getVTK().addData(m_nameScalarFields[isf], VTKFieldType::SCALAR, VTKLocation::CELL, fields[isf]);
2224 datastring.push_back(m_nameScalarFields[isf]);
2225 }
2226
2227 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
2228 int i = 0;
2229 for (const Cell &cell : mesh->getCells()) {
2230 long id = cell.getId();
2231 fieldv[ivf][i] = fieldr.vector->at(id, ivf);
2232 ++i;
2233 }
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);
2237 }
2238
2239 // Set temporary write all cells if meshPOD
2240 if (mesh == m_podkernel->getMesh())
2241 mesh->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_ALL);
2242
2243 mesh->write();
2244
2245 for (std::size_t is=0; is < datastring.size(); ++is)
2246 mesh->getVTK().removeData(datastring[is]);
2247 }
2248
2249 // Reset temporary write internal cells only if meshPOD
2250 if (mesh == m_podkernel->getMesh()){
2251#if BITPIT_ENABLE_MPI
2252 mesh->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_INTERNAL);
2253#else
2254 mesh->setVTKWriteTarget(PatchKernel::WriteTarget::WRITE_TARGET_CELLS_ALL);
2255#endif
2256 }
2257}
2258
2264void POD::readMode(std::size_t ir)
2265{
2266 if (!m_podkernel->getMesh())
2267 throw std::runtime_error ("POD: POD mesh not initialized before read a POD mode");
2268
2269 if (ir > m_nModes)
2270 throw std::runtime_error ("POD: mode index > number of POD modes during read");
2271
2272 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2273
2274 if (ir == m_nModes){
2275 // Read mean
2276 std::string filename = m_directory + m_name + ".mean"+".data";
2277 IBinaryArchive binaryReader(filename, dumpBlock);
2278 std::istream &dataStream = binaryReader.getStream();
2279
2280 // Restore solved cells (read for each POD mode)
2281 m_filter.restore(dataStream);
2282
2283 // For each mode/mean the global number of fields are read & updated again .
2284
2285 // Restore scalar fields
2286 utils::binary::read(dataStream, m_nScalarFields);
2287 m_nameScalarFields.resize(m_nScalarFields);
2288 for (std::size_t i = 0; i < m_nScalarFields; ++i)
2289 utils::binary::read(dataStream, m_nameScalarFields[i]);
2290
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);
2294 }
2295
2296 // Restore vector fields
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]);
2302
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);
2306 }
2307 binaryReader.close();
2308 }
2309 else{
2310 // Read mode
2311 m_modes[ir].clear();
2312
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();
2316
2317 // Restore solved cells (read for each POD mode)
2318 m_filter.restore(dataStream);
2319
2320 // For each mode/mean the global number of fields are read & updated again .
2321
2322 // Restore scalar fields
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]);
2327
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);
2331 }
2332
2333 // Restore vector fields
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]);
2339
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);
2343 }
2344 binaryReader.close();
2345 }
2346}
2347
2353double POD::getCellVolume(long id)
2354{
2355 return m_podkernel->getCellVolume(id);
2356}
2357
2363double POD::getRawCellVolume(long rawIndex)
2364{
2365 return m_podkernel->getRawCellVolume(rawIndex);
2366}
2367
2373void POD::dumpMode(std::size_t ir)
2374{
2375 if (ir > m_nModes)
2376 throw std::runtime_error ("POD: mode index > number of pod modes during write");
2377
2378 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2379
2380 if (ir == m_nModes) {
2381 // Dumping mean as snapshot
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);
2391 }
2392
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);
2398 }
2399 binaryWriter.close();
2400
2401 } else{
2402 // Dumping modes as snapshots
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]);
2411
2412 m_modes[ir].scalar->dump(dataStream);
2413 }
2414
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]);
2419
2420 m_modes[ir].vector->dump(dataStream);
2421 }
2422 binaryWriter.close();
2423 }
2424}
2425
2432void POD::dumpField(const std::string &name, const pod::PODField &field) const
2433{
2434 log::cout() << "Dumping snapshot ... " << std::endl;
2435
2436 int dumpBlock = (m_nProcs > 1) ? m_rank : -1;
2437
2438 // Dump field data
2439 {
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);
2444 utils::binary::write(dataStream,std::size_t(m_nScalarFields));
2445 if (m_nScalarFields) {
2446 for (std::size_t i = 0; i < m_nScalarFields; ++i) {
2447 utils::binary::write(dataStream, m_nameScalarFields[i]);
2448 }
2449 field.scalar->dump(dataStream);
2450 }
2451
2452 utils::binary::write(dataStream,std::size_t(m_nVectorFields));
2453 if (m_nVectorFields) {
2454 for (std::size_t i = 0; i < m_nVectorFields; ++i) {
2455 utils::binary::write(dataStream, m_nameVectorFields[i]);
2456 }
2457 field.vector->dump(dataStream);
2458 }
2459 binaryWriter.close();
2460 }
2461
2462 // Dumping mesh
2463 {
2464 std::string header = name;// + "mesh";
2465 OBinaryArchive binaryWriter(m_directory + name + ".mesh", ARCHIVE_VERSION, header, dumpBlock);
2466 std::ostream &dataStream = binaryWriter.getStream();
2467 field.mesh->dump(dataStream);
2468 binaryWriter.close();
2469 }
2470
2471 // Write VTK file
2472 if (m_writeMode == WriteMode::DEBUG) {
2473 write(field, name);
2474 }
2475}
2476
2483void POD::diff(pod::PODField* _a, const pod::PODMode& b)
2484{
2485 pod::PODField& a = *_a;
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");
2491 }
2492
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");
2498 }
2499
2500 if (nsf){
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;
2507 datasa++;
2508 datasb++;
2509 }
2510 }
2511 }
2512 if (nvf){
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;
2519 datava++;
2520 datavb++;
2521 }
2522 }
2523}
2524
2531void POD::sum(pod::PODField* _a, const pod::PODMode& b)
2532{
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");
2539 }
2540
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");
2546 }
2547
2548 if (nsf){
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];
2555 }
2556 }
2557 }
2558 if (nvf){
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];
2565 }
2566 }
2567 }
2568}
2569
2576std::vector<double> POD::fieldsl2norm(pod::PODField & snap)
2577{
2578 std::vector<double> norm;
2579 norm.resize(m_nFields,0.0);
2580
2581 std::unordered_set<long>::iterator it;
2582 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it){
2583 long id = *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);
2589 datas++;
2590 }
2591 }
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){
2595 norm[ifield] += dotProduct((*datav),(*datav))*snap.mesh->evalCellVolume(id);
2596 datav++;
2597 }
2598 }
2599 }
2600
2601# if BITPIT_ENABLE_MPI
2602 MPI_Allreduce(MPI_IN_PLACE, norm.data(), m_nFields, MPI_DOUBLE, MPI_SUM, m_communicator);
2603# endif
2604
2605 for (std::size_t ifield = 0; ifield < m_nFields; ++ifield)
2606 norm[ifield]=std::sqrt(norm[ifield]);
2607
2608 return norm;
2609}
2610
2617std::vector<double> POD::fieldsMax(pod::PODField & snap)
2618{
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);
2623
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));
2629 }
2630 }
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];
2635 max[i] = std::max(max[i], std::sqrt(dotProduct(fieldData, fieldData)));
2636 }
2637 }
2638 }
2639
2640# if BITPIT_ENABLE_MPI
2641 MPI_Allreduce(MPI_IN_PLACE, max.data(), m_nFields, MPI_DOUBLE, MPI_MAX, m_communicator);
2642# endif
2643
2644 return max;
2645}
2646
2647#if BITPIT_ENABLE_MPI
2653void POD::initializeCommunicator(MPI_Comm communicator)
2654{
2655 // Communication can be set just once
2656 if (isCommunicatorSet())
2657 throw std::runtime_error ("POD communicator can be set just once");
2658
2659 // The communicator has to be valid
2660 if (communicator == MPI_COMM_NULL)
2661 throw std::runtime_error ("POD communicator is not valid");
2662
2663 // Create a copy of the user-specified communicator
2664 //
2665 // No library routine should use MPI_COMM_WORLD as the communicator;
2666 // instead, a duplicate of a user-specified communicator should always
2667 // be used.
2668 MPI_Comm_dup(communicator, &m_communicator);
2669}
2670
2675MPI_Comm POD::getCommunicator() const
2676{
2677 return m_communicator;
2678}
2679
2686bool POD::isCommunicatorSet() const
2687{
2688
2689 return (getCommunicator() != MPI_COMM_NULL);
2690}
2691
2695void POD::freeCommunicator()
2696{
2697 if (!isCommunicatorSet())
2698 return;
2699
2700 int finalizedCalled;
2701 MPI_Finalized(&finalizedCalled);
2702 if (finalizedCalled)
2703 return;
2704
2705 MPI_Comm_free(&m_communicator);
2706}
2707#endif
2708
2714void POD::evalReconstructionCoeffs(pod::PODField & field)
2715{
2716 evalMinimizationMatrices();
2717
2718 // Evaluate RHS
2719 m_reconstructionCoeffs = projectField(field);
2720
2721 // Solve minimization
2722 solveMinimization(m_reconstructionCoeffs);
2723}
2724
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)
2737{
2738 std::size_t nsf = scalarIds.size();
2739 std::size_t nvf = vectorIds.size();
2740
2741 // Difference field - m_mean
2742 diff(fields, m_mean, scalarIds, podscalarIds, vectorIds, podvectorIds, &m_listActiveIDs);
2743
2744 m_reconstructionCoeffs.clear();
2745 m_reconstructionCoeffs.resize(m_nFields, std::vector<double>(m_nModes, 0.0));
2746
2747 // Evaluate minimization matrices
2748 evalMinimizationMatrices();
2749
2750 // Evaluate RHS
2751 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
2752 if (m_memoryMode == MemoryMode::MEMORY_LIGHT)
2753 readMode(ir);
2754
2755 std::unordered_set<long>::iterator it;
2756 for (it = m_listActiveIDs.begin(); it != m_listActiveIDs.end(); ++it) {
2757 long id = *it;
2758 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(id);
2759 double *datag = fields.rawData(rawIndex);
2760 if (nsf) {
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);
2766 }
2767 }
2768 if (nvf) {
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);
2775 }
2776 }
2777 }
2778 }
2779 }
2780
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);
2784#endif
2785
2786 // Solve
2787 solveMinimization(m_reconstructionCoeffs);
2788
2789 // Sum field and mean
2790 sum(fields, m_mean, scalarIds, podscalarIds, vectorIds, podvectorIds, &m_listActiveIDs);
2791}
2792
2800void POD::buildFields(const std::vector<std::vector<double>> &reconstructionCoeffs, pod::PODField &recon)
2801{
2802 // set up of the PODField members: scalar, vector
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}});
2807
2808 // Outer cycle over modes
2809 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
2810 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2811 readMode(ir);
2812 }
2813
2814 // Outer cycle over cells
2815 for (long id : m_listActiveIDs) {
2816 std::size_t rawIndex = recon.mesh->getCells().getRawIndex(id);
2817
2818 // Inner cycle over scalar fields
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];
2823 }
2824 }
2825
2826 // Inner cycle over vector fields
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];
2831 }
2832 }
2833 }
2834
2835 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2836 m_modes[ir].clear();
2837 }
2838 }
2839
2840 if (m_useMean) {
2841 sum(&recon, m_mean);
2842 }
2843}
2844
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)
2862{
2863 std::size_t nsf = scalarIds.size();
2864 std::size_t nvf = vectorIds.size();
2865 if (nsf == 0 && nvf == 0) {
2866 return;
2867 }
2868
2869 std::unordered_set<long> targetCellsStorage;
2870 if (!targetCells) {
2871 for (const Cell &cell : m_podkernel->getMesh()->getCells()) {
2872 targetCellsStorage.insert(cell.getId());
2873 }
2874 targetCells = &targetCellsStorage;
2875 }
2876
2877 // Initialization of fields
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];
2883 (*reconsi) = 0.;
2884 }
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];
2888 (*reconvi) = 0.;
2889 }
2890 }
2891 }
2892
2893 // Reconstruction of fields
2894 for (std::size_t ir = 0; ir < m_nModes; ++ir) {
2895 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2896 readMode(ir);
2897 }
2898
2899 for (const long id : *targetCells) {
2900 std::size_t rawIndex = m_podkernel->getMesh()->getCells().getRawIndex(id);
2901 double *recon = fields.rawData(rawIndex);
2902
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];
2909 }
2910 }
2911
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];
2919 }
2920 }
2921 }
2922 }
2923
2924 if (m_memoryMode == MemoryMode::MEMORY_LIGHT) {
2925 m_modes[ir].clear();
2926 }
2927 }
2928
2929 // Sum field and mean
2930 sum(fields, m_mean, scalarIds, podscalarIds, vectorIds, podvectorIds, targetCells);
2931}
2932
2939{
2940 if (!m_expert)
2941 throw std::runtime_error("POD: compute mapper can be called only in expert mode");
2942 _computeMapper(mesh);
2943}
2944
2951void POD::adaptionAlter(const std::vector<adaption::Info> & info)
2952{
2953 if (!m_expert)
2954 throw std::runtime_error("POD: update mapper can be called only in expert mode");
2955 _adaptionAlter(info);
2956}
2957
2962void POD::_computeMapper(const VolumeKernel * mesh)
2963{
2964 m_podkernel->computeMapper(mesh);
2965 m_podkernel->setMapperDirty(!m_expert);
2966}
2967
2973void POD::adaptionPrepare(const std::vector<adaption::Info> & info)
2974{
2975 m_podkernel->adaptionPrepare(info);
2976}
2977
2983void POD::_adaptionAlter(const std::vector<adaption::Info> & info)
2984{
2985 m_podkernel->adaptionAlter(info);
2986}
2987
2992void POD::adaptionCleanUp(const std::vector<adaption::Info> & info)
2993{
2994 m_podkernel->adaptionCleanUp(info);
2995}
2996
3003void POD::diff(PiercedStorage<double> &fields, const pod::PODMode &mode,
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)
3007{
3008 if (targetCells){
3009
3010 std::size_t nsf = scalarIds.size();
3011 std::size_t nvf = vectorIds.size();
3012
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);
3021 }
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];
3028 }
3029 }
3030 }
3031
3032 }
3033}
3034
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)
3045{
3046 if (targetCells){
3047
3048 std::size_t nsf = scalarIds.size();
3049 std::size_t nvf = vectorIds.size();
3050
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);
3059 }
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];
3066 }
3067 }
3068 }
3069
3070 }
3071}
3072
3079void POD::write(const pod::PODField &field, std::string file_name) const
3080{
3081 log::cout() << "Writing snapshot ... " << std::endl;
3082
3083 // Set directory to write the VTK file.
3084 m_podkernel->getMesh()->getVTK().setDirectory(m_directory);
3085
3086 // Set up empty string to store fields name
3087 std::vector<std::string > datastring;
3088 field.mesh->getVTK().setName(file_name);
3089
3090 // Set up matrix and vectors to populate the VTK file
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);
3094
3095 // Populate the scalar fields matrix
3096 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
3097 int i=0;
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);
3102 ++i;
3103 }
3104 }
3105 // Add data to the VTK file and name them using the corresponding fields name of the pod object
3106 field.mesh->getVTK().addData(m_nameScalarFields[isf], VTKFieldType::SCALAR, VTKLocation::CELL, fields[isf]);
3107 datastring.push_back(m_nameScalarFields[isf]);
3108 }
3109
3110 // Populate the vector fields matrix
3111 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
3112 int i=0;
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);
3117 ++i;
3118 }
3119 }
3120 std::string vname= m_nameVectorFields[ivf][0].substr(0,m_nameVectorFields[ivf][0].size()-2);
3121 // add data to the VTK file and name them using the corresponding fields name of the pod object
3122 field.mesh->getVTK().addData(vname, VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[ivf]);
3123 datastring.push_back(vname);
3124 }
3125
3126 // Populate the mask vector
3127 int j = 0;
3128 for (const Cell &cell : field.mesh->getCells()) {
3129 if (cell.isInterior()) {
3130 long id = cell.getId();
3131 field_mask[j] = m_filter[id];
3132 ++j;
3133 }
3134 }
3135 field.mesh->getVTK().addData("filter", VTKFieldType::SCALAR, VTKLocation::CELL, field_mask);
3136 datastring.push_back("filter");
3137
3138 // Write the VTK file
3139 field.mesh->write();
3140
3141 // Remove the data from the VTK member of the input PODField mesh
3142 for (std::size_t is=0; is < datastring.size(); ++is) {
3143 field.mesh->getVTK().removeData(datastring[is]);
3144 }
3145}
3146
3153void POD::write(int mode_index, std::string file_name)
3154{
3155 log::cout() << "Writing mode ... " << std::endl;
3156
3157 // Set up mesh
3158 VolumeKernel* mode_mesh;
3159 mode_mesh = m_podkernel->getMesh();
3160 mode_mesh->getVTK().setName(file_name);
3161
3162 // Set up empty string to store fields name
3163 std::vector<std::string> datastring;
3164
3165 // set up matrix and vectors to populate the VTK file
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);
3169
3170 // Populate the scalar fields matrix
3171 for (std::size_t isf = 0; isf < m_nScalarFields; ++isf) {
3172 int i=0;
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);
3177 ++i;
3178 }
3179 }
3180 // add data to the VTK file and name them using the corresponding fields name of the pod object
3181 mode_mesh->getVTK().addData(m_nameScalarFields[isf], VTKFieldType::SCALAR, VTKLocation::CELL, fields[isf]);
3182 datastring.push_back(m_nameScalarFields[isf]);
3183 }
3184
3185 // Populate the vector fields matrix
3186 for (std::size_t ivf = 0; ivf < m_nVectorFields; ++ivf) {
3187 int i=0;
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);
3192 ++i;
3193 }
3194 }
3195 std::string vname= m_nameVectorFields[ivf][0].substr(0,m_nameVectorFields[ivf][0].size()-2);
3196 // add data to the VTK file and name them using the corresponding fields name of the pod object
3197 mode_mesh->getVTK().addData(vname, VTKFieldType::VECTOR, VTKLocation::CELL, fieldv[ivf]);
3198 datastring.push_back(vname);
3199 }
3200
3201 // Populate the mask vector
3202 int j=0;
3203 for (const Cell &cell : mode_mesh->getCells()) {
3204 if (cell.isInterior()) {
3205 long id = cell.getId();
3206 mask[j] = m_filter[id];
3207 ++j;
3208 }
3209 }
3210 mode_mesh->getVTK().addData("mask", VTKFieldType::SCALAR, VTKLocation::CELL, mask);
3211 datastring.push_back("mask");
3212
3213 // Write the VTK file
3214 mode_mesh->write();
3215
3216 // Remove the data from the VTK member of the input PODField mesh
3217 for (std::size_t is=0; is < datastring.size(); ++is) {
3218 mode_mesh->getVTK().removeData(datastring[is]);
3219 }
3220}
3221
3222}
The Cell class defines the cells.
Definition cell.hpp:42
Input binary archive.
std::istream & getStream()
Output binary archive.
std::ostream & getStream()
The PODVolOctree is the specialized class of PODKernel for VolOctree meshes.
void setErrorThreshold(double threshold)
Definition pod.cpp:346
void computeMapper(const VolumeKernel *mesh)
Definition pod.cpp:2938
std::size_t getModeCount()
Definition pod.cpp:315
void adaptionAlter(const std::vector< adaption::Info > &info)
Definition pod.cpp:2951
void fillListActiveIDs(const PiercedStorage< bool > &bfield)
Definition pod.cpp:1315
void clear()
Definition pod.cpp:125
void run()
Definition pod.cpp:835
std::vector< std::array< std::string, 3 > > getVectorNames()
Definition pod.cpp:737
std::size_t getSnapshotCount()
Definition pod.cpp:717
void setTargetErrorFields(const std::vector< std::string > &namesf, const std::vector< std::array< std::string, 3 > > &namevf)
Definition pod.cpp:367
void evalReconstruction()
Definition pod.cpp:1448
void evalCorrelation()
Definition pod.cpp:1344
void unsetLeave1outSnapshots()
Definition pod.cpp:262
const std::string & getDirectory()
Definition pod.cpp:197
void dumpField(const std::string &name, const pod::PODField &field) const
Definition pod.cpp:2432
const std::vector< pod::PODMode > & getModes()
Definition pod.cpp:789
void setRunMode(RunMode mode)
Definition pod.cpp:541
const VolumeKernel * getMesh()
Definition pod.cpp:766
std::vector< double > fieldsMax(pod::PODField &snap)
Definition pod.cpp:2617
std::unique_ptr< PODKernel > & getKernel()
Definition pod.cpp:823
std::vector< std::string > getScalarNames()
Definition pod.cpp:727
const std::unordered_set< long int > & getListActiveIDs()
Definition pod.cpp:807
void setMeshType(MeshType type)
Definition pod.cpp:389
void setStaticMesh(bool flag)
Definition pod.cpp:482
void setSensorMask(const PiercedStorage< bool > &mask, VolumeKernel *mesh=nullptr)
Definition pod.cpp:671
WriteMode
Output Write Mode of the POD object. It defines the amount of information written by the POD object.
Definition pod.hpp:67
std::vector< double > fieldsl2norm(pod::PODField &snap)
Definition pod.cpp:2576
void setMemoryMode(MemoryMode mode)
Definition pod.cpp:504
void write(const pod::PODField &snap, std::string file_name) const
Definition pod.cpp:3079
double getErrorThreshold()
Definition pod.cpp:356
void setDirectory(const std::string &directory)
Definition pod.cpp:180
void adaptionPrepare(const std::vector< adaption::Info > &info)
Definition pod.cpp:2973
void adaptionCleanUp(const std::vector< adaption::Info > &info)
Definition pod.cpp:2992
double getEnergyLevel()
Definition pod.cpp:336
const std::string & getName()
Definition pod.cpp:170
RunMode getRunMode()
Definition pod.cpp:565
POD(MPI_Comm comm=MPI_COMM_WORLD)
Definition pod.cpp:67
void restoreDecomposition()
Definition pod.cpp:1050
void readSnapshot(const pod::SnapshotFile &snap, pod::PODField &fieldr)
Definition pod.cpp:2149
void leave1out()
Definition pod.cpp:1120
ReconstructionMode getReconstructionMode()
Definition pod.cpp:622
void dump()
Definition pod.cpp:854
std::vector< std::string > getFieldsNames()
Definition pod.cpp:748
void reconstructFields(pod::PODField &field, pod::PODField &recon)
Definition pod.cpp:1630
void evalErrorBoundingBox()
Definition pod.cpp:1502
MemoryMode getMemoryMode()
Definition pod.cpp:529
void setExpert(bool mode=true)
Definition pod.cpp:659
void setUseMean(bool flag)
Definition pod.cpp:492
void addReconstructionSnapshot(const std::string &directory, const std::string &name)
Definition pod.cpp:282
void removeLeave1outSnapshot(const std::string &directory, const std::string &name)
Definition pod.cpp:243
void evalModes()
Definition pod.cpp:2083
std::vector< std::vector< double > > projectField(pod::PODField &field)
Definition pod.cpp:1763
MeshType
Type of the Mesh used to compute the POD basis.
Definition pod.hpp:96
ReconstructionMode
Mode of Reconstruction of fields by using the POD basis.
Definition pod.hpp:77
MeshType getMeshType()
Definition pod.cpp:469
const pod::PODMode & getMean()
Definition pod.cpp:779
MemoryMode
Memory Mode of the POD object. It defines the use of the memory resources.
Definition pod.hpp:49
void setSnapshots(const std::vector< pod::SnapshotFile > &database)
Definition pod.cpp:230
void setWriteMode(WriteMode mode)
Definition pod.cpp:575
ErrorMode getErrorMode()
Definition pod.cpp:647
std::size_t getListIDInternalCount()
Definition pod.cpp:815
void setReconstructionMode(ReconstructionMode mode)
Definition pod.cpp:597
ErrorMode
Mode of Error evaluation of a reconstructed fields by the POD basis.
Definition pod.hpp:86
void dumpDecomposition()
Definition pod.cpp:903
void restore()
Definition pod.cpp:863
void setEnergyLevel(double energy)
Definition pod.cpp:326
void setModeCount(std::size_t nmodes)
Definition pod.cpp:304
void setName(const std::string &name)
Definition pod.cpp:160
WriteMode getWriteMode()
Definition pod.cpp:586
void setErrorMode(ErrorMode mode)
Definition pod.cpp:634
void evalDecomposition()
Definition pod.cpp:872
void setMesh(const std::string &directory, const std::string &name)
Definition pod.cpp:442
void evalMeanMesh()
Definition pod.cpp:1199
void evalEigen()
Definition pod.cpp:1972
std::vector< std::vector< double > > getReconstructionCoeffs()
Definition pod.cpp:799
void addSnapshot(const std::string &directory, const std::string &name)
Definition pod.cpp:208
RunMode
Run Mode of the POD object. It defines if the POD basis has to be computed or restored.
Definition pod.hpp:58
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 &)
Definition VTK.cpp:129
VTKField & addData(VTKField &&field)
Definition VTK.cpp:281
void removeData(const std::string &)
Definition VTK.cpp:330
void setName(const std::string &)
Definition VTK.cpp:119
The VolOctree defines a Octree patch.
Definition voloctree.hpp:37
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)
Definition logger.cpp:1714
Logger & info(log::Visibility defaultVisibility)
Definition logger.cpp:1856
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)
VolumeKernel * mesh
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.