37 #ifndef RBFRESCALEDSCALAR_H 38 #define RBFRESCALEDSCALAR_H 1
40 #include <lifev/core/interpolation/RBFInterpolation.hpp> 44 template <
typename mesh_Type>
97 double rbf (
double x1,
double y1,
double z1,
double x2,
double y2,
double z2,
double radius);
134 template <
typename mesh_Type>
139 template <
typename mesh_Type>
143 template <
typename mesh_Type>
153 template <
typename mesh_Type>
159 template <
typename mesh_Type>
164 M_datafile = datafile;
169 template <
typename mesh_Type>
174 this->interpolationOperator();
175 this->projectionOperator();
177 if(M_knownField->mapPtr()->commPtr()->MyPID()==0)
178 std::cout <<
"Time to assembly operators = " << TimeBuilding.diff() << std::endl;
181 this->interpolateCostantField();
184 template <
typename mesh_Type>
187 ASSERT(M_radius!=0,
"Please set the basis radius using RBFrescaledScalar<mesh_Type>::setRadius(double radius)")
189 M_neighbors.reset (
new neighbors_Type ( M_fullMeshKnown, M_localMeshKnown, M_knownField->mapPtr(), M_knownField->mapPtr()->commPtr() ) );
192 M_neighbors->setUpNeighbors();
196 M_neighbors->createPointPointNeighborsList (M_flags);
199 int LocalNodesNumber = M_GIdsKnownMesh.size();
201 std::vector<std::unordered_set<ID> > MatrixGraph (LocalNodesNumber);
202 int* ElementsPerRow =
new int[LocalNodesNumber];
203 int* GlobalID =
new int[LocalNodesNumber];
207 for (std::unordered_set<ID>::iterator it = M_GIdsKnownMesh.begin(); it != M_GIdsKnownMesh.end(); ++it)
210 MatrixGraph[k] = M_neighbors->neighborsWithinRadius (GlobalID[k], M_radius);
211 MatrixGraph[k].insert (GlobalID[k]);
212 ElementsPerRow[k] = MatrixGraph[k].size();
213 if (ElementsPerRow[k] > Max_entries)
215 Max_entries = ElementsPerRow[k];
220 M_interpolationOperatorMap.reset (
new map_Type (-1, LocalNodesNumber, GlobalID, M_knownField->mapPtr()->commPtr() ) );
221 M_interpolationOperator.reset (
new matrix_Type (*M_interpolationOperatorMap, ElementsPerRow) );
223 int* Indices =
new int[Max_entries];
224 double* Values =
new double[Max_entries];
226 for (
int i = 0 ; i < LocalNodesNumber; ++i )
229 for ( std::unordered_set<ID>::iterator it = MatrixGraph[i].begin(); it != MatrixGraph[i].end(); ++it)
232 Values[k] = rbf ( M_fullMeshKnown->point (GlobalID[i]).x(),
233 M_fullMeshKnown->point (GlobalID[i]).y(),
234 M_fullMeshKnown->point (GlobalID[i]).z(),
235 M_fullMeshKnown->point (*it).x(),
236 M_fullMeshKnown->point (*it).y(),
237 M_fullMeshKnown->point (*it).z(),
241 M_interpolationOperator->matrixPtr()->InsertGlobalValues (GlobalID[i], k, Values, Indices);
243 M_interpolationOperator->globalAssemble();
246 delete[] ElementsPerRow;
250 template <
typename mesh_Type>
256 int LocalNodesNumber = M_GIdsUnknownMesh.size();
258 std::vector<std::unordered_set<ID> > MatrixGraph (LocalNodesNumber);
259 int* ElementsPerRow =
new int[LocalNodesNumber];
260 int* GlobalID =
new int[LocalNodesNumber];
267 for (std::unordered_set<ID>::iterator it = M_GIdsUnknownMesh.begin(); it != M_GIdsUnknownMesh.end(); ++it)
271 for (
int j = 0; j < M_fullMeshKnown->numVertices(); ++j)
273 if ( M_flags[0] == -1 ||
this->isInside (M_fullMeshKnown->point (j).markerID(), M_flags) )
275 d = std::sqrt ( pow (M_fullMeshKnown->point (j).x() - M_fullMeshUnknown->point (GlobalID[k]).x(), 2)
276 + pow (M_fullMeshKnown->point (j).y() - M_fullMeshUnknown->point (GlobalID[k]).y(), 2)
277 + pow (M_fullMeshKnown->point (j).z() - M_fullMeshUnknown->point (GlobalID[k]).z(), 2) );
281 nearestPoint = M_fullMeshKnown->point (j).id();
286 MatrixGraph[k] = M_neighbors->neighborsWithinRadius (nearestPoint, M_radius);
287 MatrixGraph[k].insert (nearestPoint);
288 ElementsPerRow[k] = MatrixGraph[k].size();
289 if (ElementsPerRow[k] > Max_entries)
291 Max_entries = ElementsPerRow[k];
296 M_projectionOperatorMap.reset (
new map_Type (-1, LocalNodesNumber, GlobalID, M_unknownField->mapPtr()->commPtr() ) );
297 M_projectionOperator.reset (
new matrix_Type (*M_projectionOperatorMap, ElementsPerRow) );
299 int* Indices =
new int[Max_entries];
300 double* Values =
new double[Max_entries];
302 for (
int i = 0 ; i < LocalNodesNumber; ++i )
305 for ( std::unordered_set<ID>::iterator it = MatrixGraph[i].begin(); it != MatrixGraph[i].end(); ++it)
308 Values[k] = rbf ( M_fullMeshUnknown->point (GlobalID[i]).x(),
309 M_fullMeshUnknown->point (GlobalID[i]).y(),
310 M_fullMeshUnknown->point (GlobalID[i]).z(),
311 M_fullMeshKnown->point (*it).x(),
312 M_fullMeshKnown->point (*it).y(),
313 M_fullMeshKnown->point (*it).z(),
317 M_projectionOperator->matrixPtr()->InsertGlobalValues (GlobalID[i], k, Values, Indices);
319 M_projectionOperator->globalAssemble (M_interpolationOperatorMap, M_projectionOperatorMap);
322 delete[] ElementsPerRow;
326 template <
typename mesh_Type>
329 M_RhsF.reset (
new vector_Type (*M_interpolationOperatorMap) );
330 M_RhsOne.reset (
new vector_Type (*M_interpolationOperatorMap) );
332 M_RhsF->subset (*M_knownField, *M_interpolationOperatorMap, 0, 0);
336 template <
typename mesh_Type>
340 gamma_one.reset (
new vector_Type (*M_interpolationOperatorMap) );
346 precRawPtr->setDataFromGetPot ( M_datafile,
"prec" );
347 precPtr.reset ( precRawPtr );
349 LinearSolver solverOne;
350 solverOne.setCommunicator ( M_knownField->mapPtr()->commPtr() );
351 solverOne.setParameters ( *M_belosList );
352 solverOne.setPreconditioner ( precPtr );
355 solverOne.setRightHandSide (
M_RhsOne );
356 solverOne.solve ( gamma_one );
358 M_rbf_one.reset (
new vector_Type (*M_projectionOperatorMap) );
359 M_projectionOperator->multiply (
false, *gamma_one, *M_rbf_one);
362 template <
typename mesh_Type>
366 gamma_f.reset (
new vector_Type (*M_interpolationOperatorMap) );
372 precRawPtr->setDataFromGetPot ( M_datafile,
"prec" );
373 precPtr.reset ( precRawPtr );
375 LinearSolver solverRBF;
376 solverRBF.setCommunicator ( M_knownField->mapPtr()->commPtr() );
377 solverRBF.setParameters ( *M_belosList );
378 solverRBF.setPreconditioner ( precPtr );
381 solverRBF.setRightHandSide (
M_RhsF);
382 solverRBF.solve (gamma_f);
385 rbf_f.reset (
new vector_Type (*M_projectionOperatorMap) );
388 solution.reset (
new vector_Type (*M_projectionOperatorMap) );
390 M_projectionOperator->multiply (
false, *gamma_f, *rbf_f);
393 *solution /= *M_rbf_one;
395 M_unknownField_rbf.reset (
new vector_Type (M_unknownField->map() ) );
396 M_unknownField_rbf->subset (*rbf_f, *M_projectionOperatorMap, 0, 0);
398 M_unknownField->subset (*solution, *M_projectionOperatorMap, 0, 0);
402 template <
typename mesh_Type>
408 for ( UInt i = 0; i < LocalMesh->numVertices(); ++i )
409 if (CheckVector->blockMap().LID (
static_cast<EpetraInt_Type> (LocalMesh->point (i).id()) ) != -1)
411 GID_nodes.insert (LocalMesh->point (i).id() );
416 for ( UInt i = 0; i < LocalMesh->numVertices(); ++i )
417 if (
this->isInside (LocalMesh->point (i).markerID(), M_flags) )
418 if (CheckVector->blockMap().LID (
static_cast<EpetraInt_Type> (LocalMesh->point (i).id()) ) != -1)
420 GID_nodes.insert (LocalMesh->point (i).id() );
426 template <
typename mesh_Type>
430 for (UInt i = 0; i < flags.size(); ++i)
431 if (pointMarker == flags[i])
435 return (check > 0) ?
true :
false;
439 template <
typename mesh_Type>
440 double RBFrescaledScalar<mesh_Type>::
rbf (
double x1,
double y1,
double z1,
double x2,
double y2,
double z2,
double radius)
442 double distance = sqrt ( pow (x1 - x2, 2) + pow (y1 - y2, 2) + pow (z1 - z2, 2) );
443 return pow (1 - distance / radius, 4) * (4 * distance / radius + 1);
446 template <
typename mesh_Type>
453 template <
typename mesh_Type>
460 template <
typename mesh_Type>
467 template <
typename mesh_Type>
neighborsPtr_Type M_neighbors
vectorPtr_Type M_unknownField
std::shared_ptr< MapEpetra > mapPtr_Type
meshPtr_Type M_localMeshKnown
meshPtr_Type M_localMeshUnknown
GhostHandler< mesh_Type > neighbors_Type
void start()
Start the timer.
void solution(vectorPtr_Type &Solution)
void setupRBFData(vectorPtr_Type KnownField, vectorPtr_Type UnknownField, GetPot datafile, parameterList_Type belosList)
idContainer_Type M_GIdsUnknownMesh
idContainer_Type M_GIdsKnownMesh
static bool S_registerTetRS
std::shared_ptr< vector_Type > vectorPtr_Type
std::shared_ptr< mesh_Type > meshPtr_Type
std::shared_ptr< prec_Type > precPtr_Type
vectorPtr_Type M_unknownField_rbf
LifeV::PreconditionerIfpack prec_Type
std::vector< int > flagContainer_Type
mapPtr_Type M_interpolationOperatorMap
static const LifeV::UInt elm_nodes_num[]
mapPtr_Type M_projectionOperatorMap
matrixPtr_Type M_interpolationOperator
Epetra_Import const & importer()
Getter for the Epetra_Import.
meshPtr_Type M_fullMeshKnown
void setRadius(double radius)
void interpolateCostantField()
void interpolationOperator()
flagContainer_Type M_flags
bool isInside(ID pointMarker, flagContainer_Type Flags)
std::unordered_set< ID > idContainer_Type
RBFInterpolation< mesh_Type > * createRBFrescaledScalar()
Factory create function.
std::shared_ptr< neighbors_Type > neighborsPtr_Type
Teuchos::RCP< Teuchos::ParameterList > parameterList_Type
std::shared_ptr< matrix_Type > matrixPtr_Type
matrixPtr_Type M_projectionOperator
double rbf(double x1, double y1, double z1, double x2, double y2, double z2, double radius)
void updateRhs(vectorPtr_Type newRhs)
std::shared_ptr< basePrec_Type > basePrecPtr_Type
void solutionrbf(vectorPtr_Type &Solution_rbf)
parameterList_Type M_belosList
void identifyNodes(meshPtr_Type LocalMesh, std::unordered_set< ID > &GID_nodes, vectorPtr_Type CheckVector)
virtual ~RBFrescaledScalar()
void projectionOperator()
void stop()
Stop the timer.
void setup(meshPtr_Type fullMeshKnown, meshPtr_Type localMeshKnown, meshPtr_Type fullMeshUnknown, meshPtr_Type localMeshUnknown, flagContainer_Type flags)
meshPtr_Type M_fullMeshUnknown
LifeV::Preconditioner basePrec_Type
vectorPtr_Type M_knownField
void importFromHDF5(std::string const &fileName="ghostmap")
Import neighbor lists to an hdf5 file.
MatrixEpetra< double > matrix_Type
static bool S_registerTriRS