RBFBox.cpp
1 /*---------------------------------------------------------------------------*\
2  *
3  * mimmo
4  *
5  * Copyright (C) 2015-2021 OPTIMAD engineering Srl
6  *
7  * -------------------------------------------------------------------------
8  * License
9  * This file is part of mimmo.
10  *
11  * mimmo 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  * mimmo 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 mimmo. If not, see <http://www.gnu.org/licenses/>.
22  *
23  \ *---------------------------------------------------------------------------*/
24 
25 #include "RBFBox.hpp"
26 
27 namespace mimmo{
28 
29 
32  m_name = "mimmo.RBFBox";
33  m_origin.fill(0.0);
34  m_span.fill(0.0);
35  m_suppR = 0.0;
36 };
37 
42 RBFBox::RBFBox(const bitpit::Config::Section & rootXML){
43 
44  m_name = "mimmo.RBFBox";
45  m_origin.fill(0.0);
46  m_span.fill(0.0);
47  m_suppR = 0.0;
48 
49 
50  std::string fallback_name = "ClassNONE";
51  std::string input = rootXML.get("ClassName", fallback_name);
52  input = bitpit::utils::string::trim(input);
53  if(input == "mimmo.RBFBox"){
54  absorbSectionXML(rootXML);
55  }else{
57  };
58 }
59 
62 
66 RBFBox::RBFBox(const RBFBox & other):BaseManipulation(other){
67  m_origin = other.m_origin;
68  m_span = other.m_span;
69  m_suppR = other.m_suppR;
70 };
71 
76 void RBFBox::swap(RBFBox & x) noexcept
77 {
78  std::swap(m_origin, x.m_origin);
79  std::swap(m_span , x.m_span);
80  std::swap(m_suppR, x.m_suppR);
82 }
83 
86 void
88 
89  bool built = true;
90 
91  built = (built && createPortIn<MimmoSharedPointer<MimmoObject>, RBFBox>(this, &mimmo::RBFBox::setGeometry, M_GEOM, true));
92  built = (built && createPortIn<double, RBFBox>(this, &mimmo::RBFBox::setSupportRadius, M_VALUED));
93 
94  built = (built && createPortOut<darray3E, RBFBox>(this, &mimmo::RBFBox::getOrigin, M_POINT));
95  built = (built && createPortOut<darray3E, RBFBox>(this, &mimmo::RBFBox::getSpan, M_SPAN));
96 
97  m_arePortsBuilt = built;
98 };
99 
101 void
103  clear(); //base manipulation stuff clear
104  m_origin.fill(0.0);
105  m_span.fill(0.0);
106  m_suppR = 0.0;
107 
108 };
109 
114 darray3E
116  return(m_origin);
117 }
118 
123 darray3E
125  return(m_span);
126 }
132 void
134  bMin = m_origin - 0.5*m_span;
135  bMax = m_origin + 0.5*m_span;
136 }
137 
138 
143  if(cloud == nullptr) return;
144  if(cloud->getType() != 3) return;
146  return;
147 };
148 
152 void
154  m_suppR = std::fmax(0.0,suppR_);
155 }
156 
157 
164 void
165 RBFBox::plot(std::string directory, std::string filename,int counter, bool binary){
166 
167 #if MIMMO_ENABLE_MPI
168  if(m_rank == 0){
169 #endif
170 
171  dvecarr3E activeP(8, darray3E({{0.0,0.0,0.0}}));
172 
173  getAABB(activeP[0], activeP[6]);
174 
175  activeP[1] = activeP[0]; activeP[1][0] += m_span[0];
176  activeP[3] = activeP[0]; activeP[3][1] += m_span[1];
177  activeP[2] = activeP[6]; activeP[2][2] += -1.0*m_span[2];
178 
179  activeP[7] = activeP[6]; activeP[7][0] += -1.0*m_span[0];
180  activeP[5] = activeP[6]; activeP[5][1] += -1.0*m_span[1];
181  activeP[4] = activeP[0]; activeP[4][2] += m_span[2];
182 
183  ivector2D activeConn(1);
184  for(int i=0; i<8; ++i) activeConn[0].push_back(i);
185 
186  bitpit::VTKFormat codex = bitpit::VTKFormat::ASCII;
187  if(binary){codex=bitpit::VTKFormat::APPENDED;}
188  bitpit::VTKElementType elDM = bitpit::VTKElementType::HEXAHEDRON;
189  bitpit::VTKUnstructuredGrid vtk(directory, filename, elDM);
190  vtk.setGeomData( bitpit::VTKUnstructuredField::POINTS, activeP) ;
191  vtk.setGeomData( bitpit::VTKUnstructuredField::CONNECTIVITY, activeConn) ;
192  vtk.setDimensions(1, 8);
193  vtk.setCodex(codex);
194  if(counter>=0){vtk.setCounter(counter);}
195 
196  vtk.write();
197 #if MIMMO_ENABLE_MPI
198  }
199  MPI_Barrier(m_communicator);
200 #endif
201 };
202 
203 
208 void
210 
212  if (pc == nullptr) {
213  (*m_log)<<"Error in : "<<m_name<<" not valid point cloud linked."<<std::endl;
214  throw std::runtime_error(m_name+" not valid point cloud linked.");
215  }
216 
217  darray3E pMin, pMax;
218  pMin.fill(std::numeric_limits<double>::max());
219  pMax.fill(-1.0*std::numeric_limits<double>::max());
220 
221  int count;
222  for(const bitpit::Vertex & vert : pc->getVertices()){
223  count = 0;
224  for(double val : vert.getCoords()){
225  pMin[count] = std::min(pMin[count], (val - m_suppR*1.01));
226  pMax[count] = std::max(pMax[count], (val + m_suppR*1.01));
227  ++count;
228  }
229  }
230 
231 #if MIMMO_ENABLE_MPI
232  //reduce bbox data all over the ranks
233  MPI_Allreduce(MPI_IN_PLACE, pMin.data(), 3, MPI_DOUBLE, MPI_MIN, m_communicator);
234  MPI_Allreduce(MPI_IN_PLACE, pMax.data(), 3, MPI_DOUBLE, MPI_MAX, m_communicator);
235 #endif
236 
237  m_span = pMax - pMin;
238  m_origin = pMin + m_span * 0.5;
239 
240 };
241 
245 void
247  plot(m_outputPlot, m_name, getId(), true );
248 }
249 
255 void
256 RBFBox::absorbSectionXML(const bitpit::Config::Section & slotXML, std::string name){
257 
258  BITPIT_UNUSED(name);
259 
260  BaseManipulation::absorbSectionXML(slotXML, name);
261 
262  if(slotXML.hasOption("SupportRadius")){
263  std::string input = slotXML.get("SupportRadius");
264  input = bitpit::utils::string::trim(input);
265  double value = 0.0;
266  if(!input.empty()){
267  std::stringstream ss(input);
268  ss >> value;
269  }
270  setSupportRadius(value);
271  }
272 }
273 
279 void
280 RBFBox::flushSectionXML(bitpit::Config::Section & slotXML, std::string name){
281 
282  BITPIT_UNUSED(name);
283  BaseManipulation::flushSectionXML(slotXML, name);
284  slotXML.set("SupportRadius", std::to_string(m_suppR));
285 };
286 
287 
288 
289 }
void clearRBFBox()
Definition: RBFBox.cpp:102
#define M_GEOM
darray3E getOrigin()
Definition: RBFBox.cpp:115
void getAABB(darray3E &bMin, darray3E &bMax)
Definition: RBFBox.cpp:133
std::vector< darray3E > dvecarr3E
void buildPorts()
Definition: RBFBox.cpp:87
void setSupportRadius(double suppR_)
Definition: RBFBox.cpp:153
void setGeometry(MimmoSharedPointer< MimmoObject > cloud)
Definition: RBFBox.cpp:142
virtual void plotOptionalResults()
Definition: RBFBox.cpp:246
void plot(std::string directory, std::string filename, int counter, bool binary)
Definition: RBFBox.cpp:165
darray3E m_span
Definition: RBFBox.hpp:78
#define M_POINT
void warningXML(bitpit::Logger *log, std::string name)
BaseManipulation is the base class of any manipulation object of the library.
void execute()
Definition: RBFBox.cpp:209
Radial Basis Functions Bounding Box calculator.
Definition: RBFBox.hpp:74
virtual void absorbSectionXML(const bitpit::Config::Section &slotXML, std::string name="")
#define M_SPAN
virtual void flushSectionXML(bitpit::Config::Section &slotXML, std::string name="")
virtual void flushSectionXML(bitpit::Config::Section &slotXML, std::string name="")
Definition: RBFBox.cpp:280
void swap(RBFBox &x) noexcept
Definition: RBFBox.cpp:76
#define M_VALUED
virtual void absorbSectionXML(const bitpit::Config::Section &slotXML, std::string name="")
Definition: RBFBox.cpp:256
std::vector< ivector1D > ivector2D
void setGeometry(MimmoSharedPointer< MimmoObject > geometry)
virtual ~RBFBox()
Definition: RBFBox.cpp:61
std::array< double, 3 > darray3E
MimmoSharedPointer is a custom implementation of shared pointer.
MimmoSharedPointer< MimmoObject > getGeometry()
double m_suppR
Definition: RBFBox.hpp:79
void swap(BaseManipulation &x) noexcept
darray3E m_origin
Definition: RBFBox.hpp:77
darray3E getSpan()
Definition: RBFBox.cpp:124