Line |
Branch |
Exec |
Source |
1 |
|
|
// ______ ______ _ _ _____ ______ |
2 |
|
|
// | ____| ____| | (_)/ ____| | ____| |
3 |
|
|
// | |__ | |__ | | _| (___ ___| |__ |
4 |
|
|
// | __| | __| | | | |\___ \ / __| __| |
5 |
|
|
// | | | |____| |____| |____) | (__| |____ |
6 |
|
|
// |_| |______|______|_|_____/ \___|______| |
7 |
|
|
// Finite Elements for Life Sciences and Engineering |
8 |
|
|
// |
9 |
|
|
// License: LGL2.1 License |
10 |
|
|
// FELiScE default license: LICENSE in root folder |
11 |
|
|
// |
12 |
|
|
// Main authors: Julien Castelneau, Dominique Chapelle, Miguel Fernandez |
13 |
|
|
// |
14 |
|
|
|
15 |
|
|
// System includes |
16 |
|
|
|
17 |
|
|
// External includes |
18 |
|
|
#ifdef FELISCE_WITH_SLEPC |
19 |
|
|
#include "Core/NoThirdPartyWarning/Slepc/slepceps.hpp" |
20 |
|
|
#endif |
21 |
|
|
|
22 |
|
|
// Project includes |
23 |
|
|
#include "Solver/steklovBanner.hpp" |
24 |
|
|
#include "Solver/linearProblem.hpp" |
25 |
|
|
#include "InputOutput/io.hpp" |
26 |
|
|
#include "Core/filesystemUtil.hpp" |
27 |
|
|
|
28 |
|
|
namespace felisce { |
29 |
|
|
//! The constructor |
30 |
|
|
template<class volumeProblem> |
31 |
|
✗ |
ReducedSteklov<volumeProblem>::ReducedSteklov(PetscMatrix& mass, PetscMatrix& laplacian, MPI_Comm comm, LinearProblemReducedSteklov<volumeProblem>* pb, typename LinearProblemReducedSteklov<volumeProblem>::imgType imType, ChronoInstance::Pointer chrono): |
32 |
|
✗ |
m_comm(comm), m_mass(mass), m_laplacian(laplacian), m_lpb(pb),m_chronoRS(chrono) |
33 |
|
|
{ |
34 |
|
|
//! We need to know if we are going to extract the Dirichlet data or the Neumann data from the volume problem solution |
35 |
|
✗ |
m_imgType=imType; |
36 |
|
|
|
37 |
|
|
//! We try to load the off-line basis from file. |
38 |
|
|
|
39 |
|
|
//! The quantities: rank, nlap, numberofstekloveigenfunctions determine the folder where we hope to find the off-line basis |
40 |
|
|
//! The folder will be SteklovDataDir/NUMPROC_NBLAP_LOWRANK |
41 |
|
|
//! By default steklovDataDir = ./ |
42 |
|
|
//! This folder should be a symbolic link to another folder. |
43 |
|
✗ |
std::stringstream folder; |
44 |
|
✗ |
folder<<FelisceParam::instance().steklovDataDir |
45 |
|
✗ |
<<MpiInfo::numProc()<<"_" |
46 |
|
✗ |
<<FelisceParam::instance().nbOfLaplacianEigenFunction<<"_" |
47 |
|
✗ |
<<FelisceParam::instance().lowRank<<"/"; |
48 |
|
|
//! if the reading fails, it mains that the off-line basis is not available! |
49 |
|
✗ |
m_loaded = this->loadFromFile(folder.str()); |
50 |
|
|
} |
51 |
|
|
|
52 |
|
|
template<class volumeProblem> |
53 |
|
|
bool |
54 |
|
✗ |
ReducedSteklov<volumeProblem>::loadFromFile(std::string folder) { |
55 |
|
|
//! We try to load the following files |
56 |
|
|
//! 1. Steklov eigenVectors: |
57 |
|
|
//! eigenVectorsStreklov{0,1,2,...}.mb |
58 |
|
|
//! |
59 |
|
|
//! 2. Steklov eigenvalues: |
60 |
|
|
//! steklovEigenvalues.txt |
61 |
|
|
//! |
62 |
|
|
//! To compute these quantities you have to run the code once with some flags std::set to true (TODO) |
63 |
|
|
|
64 |
|
|
//! (1) steklov eigenVectors |
65 |
|
✗ |
for ( int i = 0; i < FelisceParam::instance().lowRank; ++i ) { |
66 |
|
✗ |
PetscVector tmp; |
67 |
|
✗ |
if (m_lpb->dofBD(/*iBD*/0).hasDofsOnBoundary() ) { |
68 |
|
✗ |
m_laplacian.getVecs(FELISCE_PETSC_NULLPTR,tmp); |
69 |
|
|
} |
70 |
|
✗ |
std::stringstream filename; |
71 |
|
✗ |
filename<<"eigenVectorsSteklov"<<i; |
72 |
|
✗ |
std::stringstream filenameBin; |
73 |
|
✗ |
filenameBin<< folder << filename.str() << ".mb"; |
74 |
|
✗ |
if ( !filesystemUtil::fileExists(filenameBin.str()) ) { |
75 |
|
✗ |
std::stringstream msg; |
76 |
|
✗ |
msg<<"Off-line basis not found in folder: "<<folder<<std::endl; |
77 |
|
✗ |
PetscPrintf(m_comm,"%s",msg.str().c_str()); |
78 |
|
✗ |
m_eVecW.clear(); |
79 |
|
✗ |
return false; |
80 |
|
|
} |
81 |
|
✗ |
if (m_lpb->dofBD(/*iBD*/0).hasDofsOnBoundary() ) { |
82 |
|
✗ |
tmp.loadFromBinaryFormat(m_comm, filename.str(), folder); |
83 |
|
|
} |
84 |
|
✗ |
m_eVecW.push_back(tmp); |
85 |
|
|
} |
86 |
|
✗ |
m_offLineBasisDim = m_eVecW.size(); |
87 |
|
✗ |
std::stringstream msg; |
88 |
|
✗ |
msg<<"Correctly loaded "<<m_offLineBasisDim<<" basis function."<<std::endl; |
89 |
|
✗ |
PetscPrintf(m_comm,"%s",msg.str().c_str()); |
90 |
|
|
|
91 |
|
|
//! (2) loading the eigenvalues |
92 |
|
✗ |
if ( ! Tools::loadVectorFromFile(m_eValues,"steklovEigenvalues.txt",folder) ) |
93 |
|
✗ |
return false; |
94 |
|
✗ |
if ( m_eValues.size() != std::size_t(FelisceParam::instance().lowRank) ) |
95 |
|
✗ |
return false; |
96 |
|
|
|
97 |
|
|
//! Finally a scaling is needed for compatibility with the rest of the code. |
98 |
|
✗ |
for ( int i(0); i < FelisceParam::instance().lowRank; ++i ) |
99 |
|
✗ |
m_eVecW[i].scale( std::sqrt(m_eValues[i]) ); |
100 |
|
|
|
101 |
|
✗ |
return true; |
102 |
|
|
|
103 |
|
|
} |
104 |
|
|
/*! \brief method that uses slepc to compute the eigenvalues/eigenvectors of the laplace-beltrami |
105 |
|
|
* operator at the interface. |
106 |
|
|
* |
107 |
|
|
* the generalized Neumann eigenproblem is solved \f[Av=\lambda Bv\f], |
108 |
|
|
* where A is neumann laplacian on the boundary, and B is the corresponding mass matrix. |
109 |
|
|
*/ |
110 |
|
|
template<class volumeProblem> |
111 |
|
|
void |
112 |
|
✗ |
ReducedSteklov<volumeProblem>::solveLaplacianEP() { |
113 |
|
|
// Necessary to compile felisce with Slepc |
114 |
|
|
#ifdef FELISCE_WITH_SLEPC |
115 |
|
|
|
116 |
|
✗ |
FEL_ASSERT(m_laplacianEVec.size()==0); |
117 |
|
|
|
118 |
|
✗ |
m_laplacianEVec.resize( FelisceParam::instance().nbOfLaplacianEigenFunction ); |
119 |
|
|
|
120 |
|
|
//! We require more eigenpairs than necessary because Slepc does not always return enough of them |
121 |
|
✗ |
int requestedNumberOfEigen = static_cast<int>( ceil( FelisceParam::instance().nbOfLaplacianEigenFunction * ( 1. + FelisceParam::instance().percentageOfExtraEig/100. ) ) ); |
122 |
|
|
|
123 |
|
|
//! We create an EPS object |
124 |
|
|
EPS laplacianEigenProblem; |
125 |
|
✗ |
EPSCreate(m_comm,&laplacianEigenProblem); |
126 |
|
✗ |
EPSSetOperators(laplacianEigenProblem,m_laplacian.toPetsc(),m_mass.toPetsc()); //<! We std::set the laplacian and the mass matrix as the two operators of the generalized problem |
127 |
|
✗ |
EPSSetProblemType(laplacianEigenProblem,EPS_GHEP); //<! We tell Slepc that the problem is GHEP: Generalized Hermitian Problem. |
128 |
|
✗ |
EPSSetDimensions(laplacianEigenProblem,requestedNumberOfEigen,PETSC_DEFAULT,PETSC_DEFAULT); |
129 |
|
✗ |
EPSSetWhichEigenpairs(laplacianEigenProblem,EPS_SMALLEST_MAGNITUDE); |
130 |
|
|
|
131 |
|
|
|
132 |
|
|
//! Necessary only when using Petsc-3.6 and only in parallel |
133 |
|
|
ST st; KSP ksp; PC pc; |
134 |
|
✗ |
EPSGetST ( laplacianEigenProblem, &st); |
135 |
|
✗ |
STGetKSP ( st, &ksp); |
136 |
|
✗ |
KSPGetPC ( ksp, &pc); |
137 |
|
✗ |
PCSetType ( pc, PCREDUNDANT); |
138 |
|
|
|
139 |
|
✗ |
EPSSetFromOptions(laplacianEigenProblem); //<! more options can be std::set in the command line. |
140 |
|
|
|
141 |
|
|
//! We solve the eigen problem. |
142 |
|
✗ |
EPSSolve(laplacianEigenProblem); |
143 |
|
✗ |
if ( FelisceParam::verbose() > 2 ) { |
144 |
|
✗ |
EPSView(laplacianEigenProblem, PETSC_VIEWER_STDOUT_(m_comm) ); |
145 |
|
|
} |
146 |
|
|
|
147 |
|
|
// We get the number of eigenpairs there were actually computed. |
148 |
|
|
PetscInt nconv; |
149 |
|
✗ |
EPSGetConverged(laplacianEigenProblem,&nconv); |
150 |
|
✗ |
PetscPrintf( m_comm," Number of converged eigenpairs: %d \n", nconv); |
151 |
|
✗ |
if ( FelisceParam::instance().nbOfLaplacianEigenFunction > nconv ) |
152 |
|
✗ |
FEL_ERROR("Number of converged eigenpairs is no enough!"); |
153 |
|
|
|
154 |
|
✗ |
std::vector<double> eigenvalues; |
155 |
|
✗ |
for ( int iEig(0); iEig < FelisceParam::instance().nbOfLaplacianEigenFunction; ++iEig ) { |
156 |
|
|
// The parallel layout of the eigen functions is taken from the laplacian matrix. |
157 |
|
✗ |
m_laplacian.getVecs(FELISCE_PETSC_NULLPTR,m_laplacianEVec[iEig]); |
158 |
|
|
|
159 |
|
|
// We extract the result of the EPS object |
160 |
|
✗ |
double eigenValue(0); |
161 |
|
✗ |
EPSGetEigenpair(laplacianEigenProblem,iEig,&eigenValue,nullptr,m_laplacianEVec[iEig].toPetsc(),nullptr); |
162 |
|
✗ |
eigenvalues.push_back(eigenValue); |
163 |
|
✗ |
if ( FelisceParam::verbose() > 1 ) { |
164 |
|
✗ |
std::stringstream master; |
165 |
|
✗ |
master<<"Laplacian eigenValues: k="<<iEig<<": "<< eigenValue << std::endl; |
166 |
|
✗ |
PetscPrintf(m_comm,"%s",master.str().c_str()); |
167 |
|
✗ |
MPI_Barrier(m_comm); |
168 |
|
|
} |
169 |
|
|
// The eigenvector related to zero is badly approximated by Slepc: we std::set it by hand |
170 |
|
✗ |
if ( eigenValue< 1e-7 ) { |
171 |
|
✗ |
std::stringstream master; |
172 |
|
✗ |
master<<"****************** Setting the "<<iEig<<" laplacian eigenValue to zero and the corresponding eigenVector to the correct std::vector."<<std::endl; |
173 |
|
✗ |
PetscPrintf(m_comm,"%s",master.str().c_str()); |
174 |
|
✗ |
MPI_Barrier(m_comm); |
175 |
|
|
|
176 |
|
✗ |
m_lpb->userSetNullSpace( m_laplacianEVec[iEig], iEig ); |
177 |
|
|
|
178 |
|
|
// Normalization |
179 |
|
✗ |
PetscVector aus; |
180 |
|
✗ |
aus.duplicateFrom( m_laplacianEVec[iEig]); |
181 |
|
✗ |
mult(m_mass, m_laplacianEVec[iEig], aus); |
182 |
|
|
double L2norm; |
183 |
|
✗ |
dot(aus,m_laplacianEVec[iEig],&L2norm); |
184 |
|
✗ |
m_laplacianEVec[iEig].scale(1./std::sqrt(L2norm) ); |
185 |
|
|
} |
186 |
|
✗ |
if ( FelisceParam::instance().exportLaplacianEigenVectors ) { |
187 |
|
✗ |
std::stringstream filename; |
188 |
|
✗ |
filename<<"laplacianEigenVector"<<iEig; |
189 |
|
✗ |
m_laplacianEVec[iEig].saveInBinaryFormat(m_comm,filename.str(),FelisceParam::instance().resultDir); |
190 |
|
|
} |
191 |
|
|
} |
192 |
|
|
// Destruction of the EPS object. |
193 |
|
✗ |
EPSDestroy(&laplacianEigenProblem); |
194 |
|
✗ |
if ( FelisceParam::instance().exportLaplacianEigenValues ) { |
195 |
|
✗ |
Tools::saveVectorInFile(eigenvalues,"laplacianEigenValues",FelisceParam::instance().resultDir); |
196 |
|
|
} |
197 |
|
|
#endif |
198 |
|
|
} |
199 |
|
|
|
200 |
|
|
|
201 |
|
|
|
202 |
|
|
/*! \brief For each eigenvector of the surface laplacian its image through the steklov operator is computed |
203 |
|
|
* |
204 |
|
|
* The result is stored in m_imgLaplacianEVec |
205 |
|
|
*/ |
206 |
|
|
template<class volumeProblem> |
207 |
|
|
void |
208 |
|
✗ |
ReducedSteklov<volumeProblem>::applySteklovOperatorOnLEV(felInt imesh) { |
209 |
|
|
|
210 |
|
✗ |
std::vector<IO::Pointer> io(1); |
211 |
|
✗ |
if ( FelisceParam::instance().exportOffLineVolumeSolution ) { |
212 |
|
✗ |
io[0] = felisce::make_shared<IO>(FelisceParam::instance().inputDirectory, FelisceParam::instance().inputFile, FelisceParam::instance().inputMesh[imesh], |
213 |
|
✗ |
FelisceParam::instance().outputMesh[imesh], FelisceParam::instance().meshDir, FelisceParam::instance().resultDir, |
214 |
|
|
"offLineSteklovSolution"); |
215 |
|
✗ |
io[0]->writeMesh(*m_lpb->mesh(imesh)); |
216 |
|
✗ |
io[0]->initializeOutput(); |
217 |
|
|
} |
218 |
|
|
|
219 |
|
✗ |
m_lpb->useSteklovDataBegin(); |
220 |
|
✗ |
FelisceParam::linearSolverVerboseOFF(); |
221 |
|
|
// Banner to improve communications on the screen |
222 |
|
✗ |
steklovBanner banner(m_lpb->m_dofBD[0/*iBD*/].numLocalDofInterface(), MpiInfo::rankProc(), m_comm, m_laplacianEVec.size() ); |
223 |
|
|
|
224 |
|
|
//The matrix 0 is cleared |
225 |
|
✗ |
m_lpb->clearMatrixRHS( FlagMatrixRHS::only_matrix ); |
226 |
|
✗ |
if ( m_lpb->numberOfMatrices() == 2 ) { |
227 |
|
|
//Also matrix 1 is cleared if there is one. |
228 |
|
✗ |
m_lpb->clearMatrix(1); |
229 |
|
|
} |
230 |
|
|
// Matrix 0 and 1 are computed and summed toghter into matrix 0. |
231 |
|
✗ |
m_lpb->assembleVolumeSystem( FlagMatrixRHS::only_matrix ); |
232 |
|
✗ |
banner.initComputation(); |
233 |
|
✗ |
for ( int iLap(0); iLap < FelisceParam::instance().nbOfLaplacianEigenFunction; ++iLap, ++banner) { |
234 |
|
✗ |
m_lpb->clearMatrixRHS(FlagMatrixRHS::only_rhs); |
235 |
|
✗ |
m_lpb->setSteklovData(m_laplacianEVec[iLap]); |
236 |
|
✗ |
m_lpb->assembleVolumeSystem(FlagMatrixRHS::only_rhs); |
237 |
|
|
|
238 |
|
✗ |
m_chronoRS->start(); |
239 |
|
✗ |
m_lpb->solve(MpiInfo::rankProc(),MpiInfo::numProc()); |
240 |
|
✗ |
m_chronoRS->stop(); |
241 |
|
|
|
242 |
|
✗ |
m_lpb->gatherSolution(); |
243 |
|
✗ |
double time=iLap; |
244 |
|
✗ |
if ( FelisceParam::instance().exportOffLineVolumeSolution ) { |
245 |
|
✗ |
m_lpb->writeSolutionFromVec(m_lpb->sequentialSolution(), MpiInfo::rankProc(), io, time, iLap, std::string("offLineVolumeSolution")); |
246 |
|
|
} |
247 |
|
✗ |
FEL_ASSERT(io[0] != nullptr) |
248 |
|
✗ |
m_imgLaplacianEVec.push_back( m_lpb->getSteklovImg() ); |
249 |
|
✗ |
if (MpiInfo::rankProc() == 0 && FelisceParam::instance().exportOffLineVolumeSolution && io[0]->typeOutput == 1 ) { //TODO and what if proc 0 has no dofs on the boundary? |
250 |
|
✗ |
io[0]->postProcess(time/*, !m_sameGeometricMeshForVariable*/); |
251 |
|
|
} |
252 |
|
|
} |
253 |
|
✗ |
banner.finalizeComputation(); |
254 |
|
|
// restoring the normal state |
255 |
|
✗ |
FelisceParam::linearSolverVerboseON(); |
256 |
|
✗ |
m_lpb->useSteklovDataEnd(); |
257 |
|
|
} |
258 |
|
|
template<class volumeProblem> |
259 |
|
|
void |
260 |
|
✗ |
ReducedSteklov<volumeProblem>::steklovEVecEstimate() { |
261 |
|
|
//! The Laplacian EigenVectors are realigned to be the eigenVectors. |
262 |
|
|
double coeff; |
263 |
|
✗ |
std::vector<PetscVector> coefficients; |
264 |
|
✗ |
this->approximateSteklovEigenVec(coefficients); |
265 |
|
|
//! The approximated eigenVectors are computed |
266 |
|
✗ |
for ( std::size_t j(0); j < coefficients.size(); ++j ) { |
267 |
|
✗ |
PetscVector currentEigenVec; |
268 |
|
✗ |
currentEigenVec.duplicateFrom( m_laplacianEVec[0]); |
269 |
|
✗ |
currentEigenVec.set(0.0); |
270 |
|
✗ |
for ( int i(0); std::size_t(i) < m_laplacianEVec.size(); ++i ) { |
271 |
|
✗ |
coefficients[j].getValues(1,&i,&coeff); |
272 |
|
✗ |
currentEigenVec.axpy(coeff,m_laplacianEVec[i]); |
273 |
|
|
} |
274 |
|
✗ |
if ( FelisceParam::instance().exportSteklovEigenVectors ) { |
275 |
|
✗ |
std::stringstream filename; |
276 |
|
✗ |
filename<<"eigenVectorsSteklov"<<j; |
277 |
|
✗ |
currentEigenVec.saveInBinaryFormat(m_comm,filename.str(),FelisceParam::instance().resultDir); |
278 |
|
|
} |
279 |
|
✗ |
currentEigenVec.scale( std::sqrt(m_eValues[j]) ); |
280 |
|
✗ |
m_eVecW.push_back(currentEigenVec); |
281 |
|
✗ |
coefficients[j].destroy(); |
282 |
|
|
} |
283 |
|
✗ |
m_offLineBasisDim = m_eVecW.size(); |
284 |
|
✗ |
m_onLineBasisDim = 0; |
285 |
|
✗ |
m_currentBasisDim = m_offLineBasisDim; |
286 |
|
|
|
287 |
|
|
//! for future loading of the off-line bais: |
288 |
|
|
|
289 |
|
|
//! rank, nlap, numberofstekloveigenfunctions determine the folder where we hope to find the off-line basis |
290 |
|
✗ |
std::stringstream folder; |
291 |
|
✗ |
folder<<FelisceParam::instance().steklovDataDir |
292 |
|
✗ |
<<MpiInfo::numProc()<<"_" |
293 |
|
✗ |
<<FelisceParam::instance().nbOfLaplacianEigenFunction<<"_" |
294 |
|
✗ |
<<FelisceParam::instance().lowRank; |
295 |
|
|
|
296 |
|
✗ |
std::stringstream createdir; |
297 |
|
✗ |
createdir<<"mkdir -p "<<FelisceParam::instance().steklovDataDir; |
298 |
|
|
|
299 |
|
✗ |
std::stringstream createlink; |
300 |
|
✗ |
std::stringstream aus; |
301 |
|
✗ |
if( FelisceParam::instance().resultDir.c_str()[0] != '/' ) |
302 |
|
✗ |
aus<<"../"<<FelisceParam::instance().resultDir; |
303 |
|
|
else |
304 |
|
✗ |
aus<<FelisceParam::instance().resultDir; |
305 |
|
✗ |
createlink<<"ln -s "<<aus.str() <<" "<<folder.str(); |
306 |
|
|
|
307 |
|
✗ |
if ( MpiInfo::rankProc() == 0 ) { |
308 |
|
|
|
309 |
|
✗ |
int ierr=system(createdir.str().c_str()); |
310 |
|
✗ |
if ( ierr>0 ) |
311 |
|
✗ |
std::cout<<"problems with: "<<createdir.str()<<std::endl; |
312 |
|
✗ |
std::cout<<createlink.str()<<std::endl; |
313 |
|
✗ |
ierr=system(createlink.str().c_str()); |
314 |
|
✗ |
if ( ierr>0 ) |
315 |
|
✗ |
std::cout<<"problems with: "<<createlink.str()<<std::endl; |
316 |
|
|
|
317 |
|
|
} |
318 |
|
|
} |
319 |
|
|
template<class volumeProblem> |
320 |
|
|
void |
321 |
|
✗ |
ReducedSteklov<volumeProblem>::assembleReducedMatrix( PetscMatrix& ReducedEigenMatrix ) { |
322 |
|
✗ |
if ( m_lpb->dofBD(/*iBD*/0).hasDofsOnBoundary() ) { |
323 |
|
|
|
324 |
|
|
//int lowRank = FelisceParam::instance().lowRank; |
325 |
|
✗ |
int nbOfLapEV = m_laplacianEVec.size(); |
326 |
|
✗ |
PetscVector tmpVec; |
327 |
|
|
|
328 |
|
|
// This matrix is supposed to be small we do not care too much about the efficiency of this code. |
329 |
|
|
//! The entries of this symmetric matrix are: (M lev_i)^T(S lev_j) |
330 |
|
✗ |
ReducedEigenMatrix.createSeqDense( PETSC_COMM_SELF, nbOfLapEV, nbOfLapEV, nullptr ); |
331 |
|
✗ |
tmpVec.duplicateFrom(m_laplacianEVec[0]); |
332 |
|
✗ |
for ( int i(0); i < nbOfLapEV; ++i ) { |
333 |
|
✗ |
std::vector<double> valuesOfCurrentRow(nbOfLapEV-i,0.0); |
334 |
|
✗ |
std::vector<int> id(nbOfLapEV-i,0); |
335 |
|
✗ |
if ( m_imgType == LinearProblemReducedSteklov<volumeProblem>::dirichlet) { |
336 |
|
|
//! M lev_i |
337 |
|
✗ |
mult(m_mass,m_laplacianEVec[i],tmpVec ); |
338 |
|
|
} else {// We do not need it in the Neumann case since it is embedded in the residual |
339 |
|
✗ |
tmpVec=m_laplacianEVec[i]; |
340 |
|
|
} |
341 |
|
|
|
342 |
|
✗ |
for ( int j(i); j < nbOfLapEV; ++j ) { |
343 |
|
✗ |
double value(0); |
344 |
|
|
//! (M lev_i)^T(S lev_j) |
345 |
|
✗ |
dot(tmpVec, m_imgLaplacianEVec[j], &value); |
346 |
|
✗ |
valuesOfCurrentRow[j-i]=value; |
347 |
|
✗ |
id[j-i]=j; |
348 |
|
|
} |
349 |
|
|
//! the matrix is symmetric: we first std::set the row and then the column |
350 |
|
✗ |
ReducedEigenMatrix.setValues( 1 , &i , nbOfLapEV-i, &id[0], &valuesOfCurrentRow[0], INSERT_VALUES); |
351 |
|
✗ |
ReducedEigenMatrix.setValues( nbOfLapEV-i-1 , &id[1] , 1, &i, &valuesOfCurrentRow[1], INSERT_VALUES); |
352 |
|
|
} |
353 |
|
✗ |
ReducedEigenMatrix.assembly(MAT_FINAL_ASSEMBLY); |
354 |
|
|
|
355 |
|
✗ |
if (FelisceParam::instance().exportReducedEigenMatrix) { |
356 |
|
✗ |
std::stringstream filename; filename<<"reducedEigenMatrix"<<MpiInfo::rankProc(); |
357 |
|
✗ |
ReducedEigenMatrix.saveInBinaryFormat(PETSC_COMM_SELF,filename.str(),FelisceParam::instance().resultDir); |
358 |
|
|
} |
359 |
|
|
|
360 |
|
|
} |
361 |
|
|
} |
362 |
|
|
/*! |
363 |
|
|
*/ |
364 |
|
|
template<class volumeProblem> |
365 |
|
|
void |
366 |
|
✗ |
ReducedSteklov<volumeProblem>::approximateSteklovEigenVec(std::vector<PetscVector> &coefficients) { |
367 |
|
|
(void) coefficients; |
368 |
|
|
#ifdef FELISCE_WITH_SLEPC |
369 |
|
✗ |
if (m_lpb->dofBD(/*iBD*/0).hasDofsOnBoundary() ) { |
370 |
|
✗ |
int lowRank = FelisceParam::instance().lowRank; |
371 |
|
|
//int nbOfLapEV = m_laplacianEVec.size(); |
372 |
|
|
|
373 |
|
|
|
374 |
|
|
// This matrix is supposed to be small we do not care too much about the efficiency of this code. |
375 |
|
|
//! The entries of this symmetric matrix are: (M lev_i)^T(S lev_j) |
376 |
|
✗ |
PetscMatrix ReducedEigenMatrix; |
377 |
|
✗ |
this->assembleReducedMatrix(ReducedEigenMatrix); |
378 |
|
|
|
379 |
|
|
//! We create an EPS object using the sub-comunicator of the procs that have at least on dof on the boundary. |
380 |
|
|
//! The matrix is supposed to be small: we do it in serial. |
381 |
|
|
EPS reducedEigenProblem; |
382 |
|
✗ |
EPSCreate(PETSC_COMM_SELF,&reducedEigenProblem); |
383 |
|
|
//! We std::set that matrix as the only operator of the eigenproblem. |
384 |
|
✗ |
EPSSetOperators(reducedEigenProblem,ReducedEigenMatrix.toPetsc(), nullptr ); |
385 |
|
|
//! We inform slepc that the problem is HEP: Hermitian Problem. |
386 |
|
✗ |
EPSSetProblemType(reducedEigenProblem,EPS_HEP); |
387 |
|
|
//! We require the eigenparis associated with the N largest eigenvalues, N=lowRank |
388 |
|
✗ |
EPSSetDimensions(reducedEigenProblem,lowRank,PETSC_DEFAULT,PETSC_DEFAULT); |
389 |
|
✗ |
if ( m_imgType == LinearProblemReducedSteklov<volumeProblem>::dirichlet) { |
390 |
|
✗ |
EPSSetWhichEigenpairs(reducedEigenProblem,EPS_LARGEST_MAGNITUDE); |
391 |
|
|
} else { |
392 |
|
✗ |
EPSSetWhichEigenpairs(reducedEigenProblem,EPS_SMALLEST_MAGNITUDE); |
393 |
|
|
} |
394 |
|
|
//EPSKrylovSchurSetRestart(reducedEigenProblem,0.5); |
395 |
|
|
//EPSSetTolerances(reducedEigenProblem,1e-8,500); |
396 |
|
✗ |
EPSSetFromOptions(reducedEigenProblem); |
397 |
|
|
|
398 |
|
|
//! We solve the eigen problem. |
399 |
|
✗ |
EPSSolve(reducedEigenProblem); |
400 |
|
✗ |
if ( FelisceParam::verbose() > 1 ) { |
401 |
|
✗ |
EPSView(reducedEigenProblem,PETSC_VIEWER_STDOUT_SELF ); |
402 |
|
|
} |
403 |
|
|
|
404 |
|
|
// We get the number of eigenpairs there were actually computed. |
405 |
|
✗ |
PetscInt nconv; EPSGetConverged(reducedEigenProblem,&nconv); |
406 |
|
✗ |
PetscPrintf(m_comm," Number of converged eigenpairs for the small eigenvalue problem: %d \n",nconv); |
407 |
|
✗ |
FEL_ASSERT( lowRank <= nconv ); |
408 |
|
✗ |
coefficients.resize(lowRank); |
409 |
|
✗ |
m_eValues.resize(lowRank); |
410 |
|
✗ |
for ( int iEig(0); iEig < lowRank; ++iEig ) { |
411 |
|
|
// The parallel layout of the eigen functions is taken from the operator of the problem. |
412 |
|
✗ |
ReducedEigenMatrix.getVecs(FELISCE_PETSC_NULLPTR,coefficients[iEig]); |
413 |
|
|
|
414 |
|
|
// We extract the result of the EPS object |
415 |
|
✗ |
EPSGetEigenpair(reducedEigenProblem,iEig,&m_eValues[iEig],nullptr,coefficients[iEig].toPetsc(),nullptr); |
416 |
|
✗ |
if ( FelisceParam::verbose() > 1) { |
417 |
|
✗ |
std::stringstream master; |
418 |
|
✗ |
master<<"Eigen Values: "<<iEig<<" , "<< m_eValues[iEig] << std::endl; |
419 |
|
✗ |
PetscPrintf(m_comm,"%s",master.str().c_str()); |
420 |
|
✗ |
MPI_Barrier(m_comm); |
421 |
|
|
} |
422 |
|
✗ |
if (FelisceParam::instance().exportFTOfSteklovEigenVectors && MpiInfo::rankProc() == 0 ) { |
423 |
|
✗ |
std::stringstream filename; filename<<"fourierTransformOfSteklovEigenVectors"<<iEig; |
424 |
|
✗ |
coefficients[iEig].saveInBinaryFormat(PETSC_COMM_SELF,filename.str(),FelisceParam::instance().resultDir); |
425 |
|
|
} |
426 |
|
|
} |
427 |
|
|
// destruction of the EPS object. |
428 |
|
✗ |
EPSDestroy(&reducedEigenProblem); |
429 |
|
✗ |
ReducedEigenMatrix.destroy(); |
430 |
|
✗ |
if (FelisceParam::instance().exportSteklovEigenValues) { |
431 |
|
✗ |
Tools::saveVectorInFile(m_eValues,"steklovEigenvalues.txt",FelisceParam::instance().resultDir); |
432 |
|
|
} |
433 |
|
|
} |
434 |
|
|
#else |
435 |
|
|
std::cerr<<"You have to compile this program with the flag FELISCE_WITH_SLEPC"<< std::endl; |
436 |
|
|
#endif |
437 |
|
|
} |
438 |
|
|
|
439 |
|
|
|
440 |
|
|
/*! \brief The action of the low rank steklov on an a given input std::vector |
441 |
|
|
* |
442 |
|
|
* \param[in] in, the given std::vector |
443 |
|
|
* \param[out] out, the outcome: it is supposed to be already allocated. |
444 |
|
|
* |
445 |
|
|
* We use the formula \f[ S = V L V^T M \f], where V is matrix with the eigenvectors, L of the eigenvalues and M is the mass. |
446 |
|
|
* If the tryAcceleration parameter is on: |
447 |
|
|
* First we check if the input is well represented on the current basis, if not it is added to the basis after having solved the volume problem |
448 |
|
|
*/ |
449 |
|
|
template<class volumeProblem> |
450 |
|
|
void |
451 |
|
✗ |
ReducedSteklov<volumeProblem>::applyLowRankSteklov(PetscVector& in, PetscVector& out) { |
452 |
|
✗ |
m_chronoRS->start(); |
453 |
|
|
// Initialization phase, out reset to zero since we comulate on it in the different for loops |
454 |
|
✗ |
out.set(0.0); |
455 |
|
|
|
456 |
|
|
// Auxiliary double |
457 |
|
✗ |
double partialResult(0); |
458 |
|
|
|
459 |
|
✗ |
PetscVector tmpVec,ausForNeu; |
460 |
|
✗ |
tmpVec.duplicateFrom(in); |
461 |
|
✗ |
ausForNeu.duplicateFrom(in); |
462 |
|
|
// The quantity M*in is stored since it will used several times |
463 |
|
✗ |
mult(m_mass,in,tmpVec ); |
464 |
|
|
|
465 |
|
|
|
466 |
|
✗ |
if ( ! FelisceParam::instance().tryAcceleration ) { |
467 |
|
✗ |
for ( std::size_t i(0); i < m_offLineBasisDim; ++i ) { |
468 |
|
✗ |
if ( ! FelisceParam::instance().onlyConstResp ) { |
469 |
|
|
//! \f[ partialResult = q_i \cdot M in \f] |
470 |
|
✗ |
dot( m_eVecW[i], tmpVec, &partialResult); |
471 |
|
✗ |
if ( m_imgType == LinearProblemReducedSteklov<volumeProblem>::neumann) { |
472 |
|
|
// We pre-multiply the result by -M because we want to pass the equivalent of a residual to the coupled problem |
473 |
|
✗ |
partialResult *= -1; |
474 |
|
✗ |
mult(m_mass, m_eVecW[i], ausForNeu); |
475 |
|
✗ |
out.axpy( partialResult, ausForNeu); //out = - pr M q (we want the residual not the traction ) //DEBUG |
476 |
|
|
} else { |
477 |
|
|
//! \f[ out += partialResult q_i \f] |
478 |
|
✗ |
out.axpy( partialResult, m_eVecW[i]); |
479 |
|
|
} |
480 |
|
|
} |
481 |
|
|
} |
482 |
|
|
} else { |
483 |
|
✗ |
if ( m_imgType == LinearProblemReducedSteklov<volumeProblem>::neumann) { |
484 |
|
✗ |
std::cout<<"WARNING: this part has not been checked for dirichlet 2 neumann problem"<<std::endl; |
485 |
|
✗ |
std::cout<<__FILE__<<":"<<__LINE__<<std::endl; |
486 |
|
|
} |
487 |
|
✗ |
std::size_t curOnLineDim = m_onLineBasisDim; |
488 |
|
|
//! 1) COMPUTATION OF THE TRUE L2-NORM OF THE INPUT |
489 |
|
|
// L2 norm of the input std::vector |
490 |
|
✗ |
double L2Norm(0); |
491 |
|
✗ |
dot(in,tmpVec, &L2Norm); |
492 |
|
✗ |
L2Norm=std::sqrt(L2Norm); |
493 |
|
|
|
494 |
|
|
//! 2) COMPUTATION OF THE FOURIER COEFFICIENTS |
495 |
|
|
//! 2bis) COMPUTATION OF THE L2-NORM OF THE INPUT'S PROJECTION |
496 |
|
|
// now we compute the norm of the projection onto the basis |
497 |
|
✗ |
double projNorm(0); |
498 |
|
|
// Vector that'll contains all the scalar product w.r.t. the vectors of the basis |
499 |
|
✗ |
std::vector<double> partialResults( m_currentBasisDim ); |
500 |
|
|
// loop on the offline basis |
501 |
|
✗ |
for ( std::size_t i(0); i < m_offLineBasisDim; ++i ) { |
502 |
|
|
// be careful: the eigenVector contains the square root of the eigenValue |
503 |
|
✗ |
dot(m_eVecW[i], tmpVec, &partialResult); |
504 |
|
✗ |
partialResults[i]=partialResult; |
505 |
|
✗ |
projNorm += std::pow( partialResult / std::sqrt( m_eValues[i] ), 2 ); |
506 |
|
|
} |
507 |
|
✗ |
for ( std::size_t i(0); i< m_onLineBasisDim; ++i ) { |
508 |
|
✗ |
dot(m_onlineInputs[i],tmpVec, &partialResult); |
509 |
|
✗ |
partialResults[ i + m_offLineBasisDim ] = partialResult; |
510 |
|
✗ |
projNorm += std::pow( partialResult, 2 ); |
511 |
|
|
} |
512 |
|
✗ |
projNorm=std::sqrt(projNorm); |
513 |
|
|
|
514 |
|
|
//! 3) COMPUTATION OF THE RATIO: ||P(X_in)||/|| X_in || |
515 |
|
✗ |
double ratio(projNorm/L2Norm); |
516 |
|
✗ |
if( MpiInfo::rankProc() == 0 ) |
517 |
|
✗ |
std::cout<<"percentage of input vectors that is taken into account by the approximation: "<<ratio*100<<std::endl; |
518 |
|
|
|
519 |
|
✗ |
if ( ratio < FelisceParam::instance().projOnlineRatio )//0.995 |
520 |
|
|
{ |
521 |
|
✗ |
if( MpiInfo::rankProc() == 0 ) |
522 |
|
✗ |
std::cout<<"Too complicated input, the volume problem will be solved!"<<std::endl; |
523 |
|
|
//! 4) CONSTRUCTION OF A NEW BASIS VECTOR |
524 |
|
|
// Reconstruct the projection of in, compute in - P(in) |
525 |
|
✗ |
PetscVector newBasisVector; |
526 |
|
|
// first we take In and then we remove its projections on the current basis |
527 |
|
✗ |
newBasisVector.duplicateFrom(in); |
528 |
|
✗ |
newBasisVector.copyFrom(in); |
529 |
|
✗ |
for ( std::size_t i(0); i< m_offLineBasisDim; ++i ) { |
530 |
|
✗ |
newBasisVector.axpy( - partialResults[i] / m_eValues[i], m_eVecW[i]); |
531 |
|
|
} |
532 |
|
✗ |
for ( std::size_t i(0); i<m_onLineBasisDim; ++i ) { |
533 |
|
✗ |
newBasisVector.axpy( - partialResults[ i + m_eVecW.size() ], m_onlineInputs[i] ); |
534 |
|
|
} |
535 |
|
|
// renormalize it, with respect to Mass (we know the norm by difference) |
536 |
|
✗ |
newBasisVector.scale( 1./( std::sqrt(L2Norm*L2Norm - projNorm*projNorm) ) ); |
537 |
|
|
|
538 |
|
|
// DEBUG |
539 |
|
|
#ifdef FELISCE_DEBUG |
540 |
|
|
double check; |
541 |
|
✗ |
PetscVector tmpVec2; |
542 |
|
|
//Computing L2-norm |
543 |
|
✗ |
tmpVec2.duplicateFrom(newBasisVector); |
544 |
|
✗ |
mult(m_mass,newBasisVector,tmpVec2 ); |
545 |
|
✗ |
dot(tmpVec2, newBasisVector, &check); |
546 |
|
✗ |
tmpVec2.destroy(); |
547 |
|
|
//checking that it has been normalized |
548 |
|
✗ |
if ( MpiInfo::rankProc() == 0 ) |
549 |
|
✗ |
std::cout<<"added a new input Vector, norm should be one and it is: "<< check << std::endl; |
550 |
|
✗ |
FEL_ASSERT( Tools::almostEqual(check,1.,1.e-12) ); |
551 |
|
|
#endif |
552 |
|
|
//! 5) ADDITION OF THE NEW VECTOR TO THE ONLINE BASIS |
553 |
|
✗ |
m_onlineInputs.push_back(newBasisVector); |
554 |
|
✗ |
++m_onLineBasisDim; |
555 |
|
✗ |
++m_currentBasisDim; |
556 |
|
|
|
557 |
|
|
//! 6) COMPUTATION OF THE IMAGE OF THE NEW VECTOR |
558 |
|
✗ |
m_lpb->useSteklovDataBegin(); |
559 |
|
✗ |
m_lpb->setSteklovData( newBasisVector ); |
560 |
|
✗ |
m_lpb->clearMatrixRHS(FlagMatrixRHS::matrix_and_rhs); |
561 |
|
✗ |
m_lpb->assembleVolumeSystem(FlagMatrixRHS::matrix_and_rhs); |
562 |
|
✗ |
m_lpb->solve(MpiInfo::rankProc(),MpiInfo::numProc()); |
563 |
|
✗ |
out = m_lpb->getSteklovImg(); |
564 |
|
✗ |
m_lpb->useSteklovDataEnd(); |
565 |
|
|
//! 7) SAVING OF THE IMAGE OF THE NEW VECTOR |
566 |
|
|
// save in a list along with its image. |
567 |
|
✗ |
PetscVector imageOfNewBasisVec; |
568 |
|
✗ |
imageOfNewBasisVec.duplicateFrom(out); |
569 |
|
✗ |
imageOfNewBasisVec.copyFrom(out); |
570 |
|
✗ |
m_onlineOutputs.push_back(imageOfNewBasisVec); |
571 |
|
|
//! 8) COMPUTATION OF THE CURRENT IMAGE ON THE NEW BASIS |
572 |
|
|
// the output is rescaled with the correct norm |
573 |
|
✗ |
out.scale(std::sqrt(L2Norm*L2Norm - projNorm*projNorm)); |
574 |
|
|
} |
575 |
|
|
//! 9) COMPUTATION OF THE CURRENT IMAGE WITH THE REST OF THE BASIS |
576 |
|
✗ |
for ( std::size_t i(0); i < m_offLineBasisDim; ++i ) { |
577 |
|
✗ |
out.axpy( partialResults[i], m_eVecW[i]); |
578 |
|
|
} |
579 |
|
✗ |
for ( std::size_t i(0); i < curOnLineDim; ++i ) { |
580 |
|
✗ |
out.axpy( partialResults[ i + m_offLineBasisDim ], m_onlineOutputs[i]); |
581 |
|
|
} |
582 |
|
|
} |
583 |
|
✗ |
if ( m_constantResponse.size() == 1) { |
584 |
|
|
//TODO: case with multiple constant responses |
585 |
|
✗ |
double coeff = m_lpb->userInletPressure() * ( ( m_imgType == LinearProblemReducedSteklov<volumeProblem>::neumann ) ? -1 : 1 ); |
586 |
|
✗ |
out.axpy( coeff,m_constantResponse[0]); |
587 |
|
|
} |
588 |
|
✗ |
m_chronoRS->stop(); |
589 |
|
|
} |
590 |
|
|
|
591 |
|
|
template<class volumeProblem> |
592 |
|
|
void |
593 |
|
✗ |
ReducedSteklov<volumeProblem>::exportAllEig( felInt imesh ) { |
594 |
|
|
|
595 |
|
✗ |
std::vector<IO::Pointer> ioSteklov = {felisce::make_shared<IO>(FelisceParam::instance().inputDirectory, FelisceParam::instance().inputFile, FelisceParam::instance().inputMesh[imesh], |
596 |
|
✗ |
FelisceParam::instance().outputMesh[imesh], FelisceParam::instance().meshDir, FelisceParam::instance().resultDir, |
597 |
|
|
"steklovEigen")}; |
598 |
|
|
|
599 |
|
✗ |
std::vector<IO::Pointer> ioLaplacian = {felisce::make_shared<IO>(FelisceParam::instance().inputDirectory, FelisceParam::instance().inputFile, FelisceParam::instance().inputMesh[imesh], |
600 |
|
✗ |
FelisceParam::instance().outputMesh[imesh], FelisceParam::instance().meshDir, FelisceParam::instance().resultDir, |
601 |
|
|
"laplacianEigen")}; |
602 |
|
|
|
603 |
|
✗ |
ioSteklov[0] ->writeMesh(*m_lpb->m_mesh[0]); |
604 |
|
✗ |
ioLaplacian[0] ->writeMesh(*m_lpb->m_mesh[0]); |
605 |
|
|
|
606 |
|
✗ |
ioSteklov[0] ->initializeOutput(); |
607 |
|
✗ |
ioLaplacian[0] ->initializeOutput(); |
608 |
|
|
|
609 |
|
✗ |
PetscVector ausPara,ausSeq; |
610 |
|
✗ |
ausPara.duplicateFrom(m_lpb->solution()); |
611 |
|
✗ |
ausSeq.duplicateFrom(m_lpb->sequentialSolution()); |
612 |
|
|
|
613 |
|
✗ |
int nEVECWl = m_eVecW.size(); |
614 |
|
✗ |
int nEVECW(0); |
615 |
|
✗ |
MPI_Allreduce( &nEVECWl,&nEVECW,1,MPI_INT,MPI_MAX,MpiInfo::petscComm()); |
616 |
|
|
|
617 |
|
✗ |
if ( m_eValues.size() < (std::size_t)nEVECW ) |
618 |
|
✗ |
m_eValues.resize(nEVECW); |
619 |
|
|
|
620 |
|
✗ |
for ( std::size_t k(0); k<(std::size_t)nEVECW; ++k) { |
621 |
|
✗ |
m_lpb->dofBD(/*iBD*/0).extendOnVolume( ausPara, m_eVecW[k]); |
622 |
|
✗ |
m_lpb->gatherVector(ausPara,ausSeq); |
623 |
|
✗ |
ausSeq.scale( 1./std::sqrt( m_eValues[k] ) ); |
624 |
|
✗ |
double time=k; |
625 |
|
✗ |
m_lpb->writeSolutionFromVec(ausSeq, MpiInfo::rankProc(), ioSteklov, time, k, std::string("steklovEigenVec")); |
626 |
|
✗ |
if (MpiInfo::rankProc() == 0) { |
627 |
|
✗ |
if ( ioSteklov[0]->typeOutput == 1 ) // 1 : ensight{ |
628 |
|
✗ |
ioSteklov[0]->postProcess(time/*, !m_sameGeometricMeshForVariable*/); |
629 |
|
|
} |
630 |
|
|
} |
631 |
|
|
|
632 |
|
|
|
633 |
|
✗ |
int nLAPl = m_laplacianEVec.size(); |
634 |
|
✗ |
int nLAP(0); |
635 |
|
✗ |
MPI_Allreduce( &nLAPl,&nLAP,1,MPI_INT,MPI_MAX,MpiInfo::petscComm()); |
636 |
|
|
|
637 |
|
✗ |
for ( std::size_t k(0); k<(std::size_t)nLAP; ++k) { |
638 |
|
✗ |
m_lpb->dofBD(/*iBD*/0).extendOnVolume( ausPara, m_laplacianEVec[k]); |
639 |
|
✗ |
m_lpb->gatherVector(ausPara,ausSeq); |
640 |
|
✗ |
double time=k; |
641 |
|
✗ |
m_lpb->writeSolutionFromVec(ausSeq, MpiInfo::rankProc(), ioLaplacian, time, k, std::string("laplacianEigenVec")); |
642 |
|
✗ |
if (MpiInfo::rankProc() == 0) { |
643 |
|
✗ |
if ( ioLaplacian[0]->typeOutput == 1 ) // 1 : ensight{ |
644 |
|
✗ |
ioLaplacian[0]->postProcess(time/*, !m_sameGeometricMeshForVariable*/); |
645 |
|
|
} |
646 |
|
|
} |
647 |
|
|
} |
648 |
|
|
} |
649 |
|
|
|
650 |
|
|
|