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 |
|
|
#include <stack> |
17 |
|
|
|
18 |
|
|
// External includes |
19 |
|
|
|
20 |
|
|
// Project includes |
21 |
|
|
#include "linearProblemNSContinuation.hpp" |
22 |
|
|
#include "Core/felisceTransient.hpp" |
23 |
|
|
#include "FiniteElement/elementVector.hpp" |
24 |
|
|
#include "FiniteElement/elementMatrix.hpp" |
25 |
|
|
#include "InputOutput/io.hpp" |
26 |
|
|
#include "DegreeOfFreedom/boundaryCondition.hpp" |
27 |
|
|
|
28 |
|
|
/*! |
29 |
|
|
\file linearProblemNSContinuation.cpp |
30 |
|
|
\date 22/02/2022 |
31 |
|
|
\brief Continuation method for FSI, monolithic |
32 |
|
|
*/ |
33 |
|
|
|
34 |
|
|
namespace felisce { |
35 |
|
|
|
36 |
|
✗ |
LinearProblemNSContinuation::LinearProblemNSContinuation(): |
37 |
|
✗ |
LinearProblem() |
38 |
|
|
{ |
39 |
|
✗ |
m_useDataPreviousVelocity = true ; |
40 |
|
✗ |
m_previousVelocity = &m_oldVelocity; |
41 |
|
✗ |
m_aoPreviousVelocity = &m_ao; |
42 |
|
✗ |
m_dofPreviousVelocity = &m_dof; |
43 |
|
|
} |
44 |
|
|
|
45 |
|
✗ |
void LinearProblemNSContinuation::setPreviousVelocity(PetscVector & previousVelocityCandidate, AO & aoPreviousVelocity, Dof & dofPreviousVelocity) |
46 |
|
|
{ |
47 |
|
✗ |
m_previousVelocity = &previousVelocityCandidate; |
48 |
|
✗ |
m_aoPreviousVelocity = &aoPreviousVelocity; |
49 |
|
✗ |
m_dofPreviousVelocity = &dofPreviousVelocity; |
50 |
|
|
|
51 |
|
✗ |
m_useDataPreviousVelocity = false ; |
52 |
|
|
} |
53 |
|
|
|
54 |
|
✗ |
void LinearProblemNSContinuation::readVelocityData(IO& io,double iteration) { |
55 |
|
|
// Read velocity variable from the .case file |
56 |
|
|
// store previous velocity from data |
57 |
|
✗ |
if(m_useDataPreviousVelocity) // validation |
58 |
|
|
{ |
59 |
|
✗ |
if(iteration==0 || (m_velocityData.isNull()) ) |
60 |
|
|
{ |
61 |
|
✗ |
m_oldVelocity.duplicateFrom(this->sequentialSolution()); |
62 |
|
✗ |
m_oldVelocity.set(0.); |
63 |
|
|
}else{ |
64 |
|
✗ |
m_oldVelocity.copyFrom(m_velocityData); |
65 |
|
|
} |
66 |
|
|
} |
67 |
|
|
|
68 |
|
✗ |
if(iteration==0 || (m_velocityData.isNull()) ) |
69 |
|
✗ |
m_velocityData.duplicateFrom(this->sequentialSolution()); |
70 |
|
|
|
71 |
|
✗ |
std::vector<double> dataVelocity; |
72 |
|
✗ |
dataVelocity.resize(m_mesh[m_currentMesh]->numPoints()*3, 0.); |
73 |
|
✗ |
io.readVariable(0,iteration, dataVelocity.data(), dataVelocity.size()); |
74 |
|
|
|
75 |
|
✗ |
int dim = this->dimension(); |
76 |
|
|
|
77 |
|
✗ |
m_vectorVelocityData.resize(m_mesh[m_currentMesh]->numPoints()*dim, 0.); |
78 |
|
|
|
79 |
|
|
|
80 |
|
✗ |
for (int i=0; i < m_mesh[m_currentMesh]->numPoints(); ++i) { |
81 |
|
✗ |
for(int icomp=0; icomp<dim; icomp++) { |
82 |
|
✗ |
m_vectorVelocityData[dim*i+icomp]=dataVelocity[3*i+icomp]; |
83 |
|
|
} |
84 |
|
|
} |
85 |
|
|
|
86 |
|
✗ |
std::vector<PetscInt> ao_potData( m_vectorVelocityData.size()); |
87 |
|
✗ |
for (size_t i=0; i < m_vectorVelocityData.size(); ++i) { |
88 |
|
✗ |
ao_potData[i]=i; |
89 |
|
|
} |
90 |
|
|
|
91 |
|
✗ |
AOApplicationToPetsc(this->ao(),m_vectorVelocityData.size(),ao_potData.data()); |
92 |
|
✗ |
VecSetValues(m_velocityData.toPetsc(),m_vectorVelocityData.size(),ao_potData.data(), m_vectorVelocityData.data(), INSERT_VALUES); |
93 |
|
|
|
94 |
|
|
//m_velocityData.view(PETSC_VIEWER_STDOUT_SELF); |
95 |
|
|
} |
96 |
|
|
|
97 |
|
|
|
98 |
|
✗ |
void LinearProblemNSContinuation::initialize(std::vector<GeometricMeshRegion::Pointer>& mesh, FelisceTransient::Pointer fstransient, MPI_Comm& comm, bool doUseSNES) { |
99 |
|
✗ |
LinearProblem::initialize(mesh,comm, doUseSNES); |
100 |
|
✗ |
this->m_fstransient = fstransient; |
101 |
|
✗ |
std::vector<PhysicalVariable> listVariable; |
102 |
|
✗ |
std::vector<std::size_t> listNumComp; |
103 |
|
|
|
104 |
|
|
// primal variables |
105 |
|
✗ |
listVariable.push_back(velocity); |
106 |
|
✗ |
listNumComp.push_back(this->dimension()); |
107 |
|
✗ |
listVariable.push_back(pressure); |
108 |
|
✗ |
listNumComp.push_back(1); |
109 |
|
|
|
110 |
|
|
// dual variables with dummy names |
111 |
|
✗ |
listVariable.push_back(velocityDiv); |
112 |
|
✗ |
listNumComp.push_back(this->dimension()); |
113 |
|
✗ |
listVariable.push_back(potThorax); |
114 |
|
✗ |
listNumComp.push_back(1); |
115 |
|
|
|
116 |
|
✗ |
if(FelisceParam::instance(this->instanceIndex()).useALEformulation) |
117 |
|
|
{ |
118 |
|
✗ |
listVariable.push_back(displacement); |
119 |
|
✗ |
listNumComp.push_back(this->dimension()); |
120 |
|
|
} |
121 |
|
|
|
122 |
|
|
|
123 |
|
|
//define unknown of the linear system. |
124 |
|
✗ |
m_listUnknown.push_back(velocity); |
125 |
|
✗ |
m_listUnknown.push_back(pressure); |
126 |
|
✗ |
m_listUnknown.push_back(velocityDiv); |
127 |
|
✗ |
m_listUnknown.push_back(potThorax); |
128 |
|
|
|
129 |
|
✗ |
userAddOtherUnknowns(listVariable,listNumComp); |
130 |
|
✗ |
userAddOtherVariables(listVariable,listNumComp); |
131 |
|
✗ |
definePhysicalVariable(listVariable,listNumComp); |
132 |
|
|
|
133 |
|
✗ |
m_iVelocity = m_listVariable.getVariableIdList(velocity); |
134 |
|
✗ |
m_iPressure = m_listVariable.getVariableIdList(pressure); |
135 |
|
✗ |
m_iVelocityDiv = m_listVariable.getVariableIdList(velocityDiv); |
136 |
|
✗ |
m_iPotThorax = m_listVariable.getVariableIdList(potThorax); |
137 |
|
|
|
138 |
|
✗ |
m_viscosity = FelisceParam::instance(this->instanceIndex()).viscosity; |
139 |
|
✗ |
m_useSymmetricStress = FelisceParam::instance(this->instanceIndex()).useSymmetricStress; |
140 |
|
|
|
141 |
|
|
|
142 |
|
|
} |
143 |
|
|
|
144 |
|
✗ |
void LinearProblemNSContinuation::initPerElementType(ElementType eltType, FlagMatrixRHS flagMatrixRHS) { |
145 |
|
|
IGNORE_UNUSED_ARGUMENT(eltType); |
146 |
|
|
IGNORE_UNUSED_ARGUMENT(flagMatrixRHS); |
147 |
|
✗ |
m_feVel = m_listCurrentFiniteElement[m_iVelocity]; |
148 |
|
✗ |
m_fePres = m_listCurrentFiniteElement[m_iPressure]; |
149 |
|
|
|
150 |
|
|
// for the convective term + ALE |
151 |
|
✗ |
m_elemFieldAdv.initialize(DOF_FIELD, *m_feVel, this->dimension()); |
152 |
|
|
|
153 |
|
|
// define elemField to store the given velocity measurements and the source terms |
154 |
|
✗ |
m_dataTerm.initialize(DOF_FIELD, *m_feVel, this->dimension()); |
155 |
|
✗ |
m_sourceTerm.initialize(DOF_FIELD, *m_feVel, this->dimension()); |
156 |
|
|
|
157 |
|
|
//========= |
158 |
|
|
// ALE |
159 |
|
|
//========= |
160 |
|
✗ |
if(FelisceParam::instance(this->instanceIndex()).useALEformulation) { |
161 |
|
✗ |
m_iDisplacement = m_listVariable.getVariableIdList(displacement); |
162 |
|
|
|
163 |
|
✗ |
m_elemFieldVelMesh.initialize(DOF_FIELD, *m_feVel, this->dimension()); |
164 |
|
✗ |
m_beta.initialize(DOF_FIELD, *m_feVel, this->dimension()); |
165 |
|
|
|
166 |
|
✗ |
numDofExtensionProblem = m_externalDof[0]->numDof(); |
167 |
|
|
|
168 |
|
✗ |
m_petscToGlobal1.resize(numDofExtensionProblem); |
169 |
|
✗ |
m_petscToGlobal2.resize(numDofExtensionProblem); |
170 |
|
✗ |
m_auxvec.resize(numDofExtensionProblem); |
171 |
|
|
|
172 |
|
✗ |
for (int i=0; i<numDofExtensionProblem; i++) { |
173 |
|
✗ |
m_petscToGlobal1[i]=i; |
174 |
|
✗ |
m_petscToGlobal2[i]=i; |
175 |
|
|
} |
176 |
|
✗ |
AOApplicationToPetsc(m_externalAO[0],numDofExtensionProblem,m_petscToGlobal1.data()); |
177 |
|
✗ |
AOApplicationToPetsc(m_ao,numDofExtensionProblem,m_petscToGlobal2.data()); |
178 |
|
|
} |
179 |
|
|
} |
180 |
|
|
|
181 |
|
|
|
182 |
|
✗ |
void LinearProblemNSContinuation::computeElementArray(const std::vector<Point*>& elemPoint, const std::vector<felInt>& /*elemIdPoint*/, felInt& iel, FlagMatrixRHS flagMatrixRHS) { |
183 |
|
✗ |
m_feVel->updateFirstDeriv(0, elemPoint); |
184 |
|
✗ |
m_fePres->updateFirstDeriv(0, elemPoint); |
185 |
|
|
|
186 |
|
✗ |
double hK = m_feVel->diameter(); |
187 |
|
|
|
188 |
|
|
// stabilization parameters |
189 |
|
✗ |
double gamma_M = FelisceParam::instance(this->instanceIndex()).referenceValue2; // data fidelity coefficient |
190 |
|
✗ |
double gamma_p = hK*hK*FelisceParam::instance(this->instanceIndex()).referenceValue3; // Brezzi-Pitkaranta coefficient |
191 |
|
✗ |
double gamma_p_star = FelisceParam::instance(this->instanceIndex()).referenceValue4; // dual pressure coefficient |
192 |
|
✗ |
double gamma_div = 0.1; |
193 |
|
✗ |
double gamma_u_star = 0.1; |
194 |
|
|
|
195 |
|
|
// Luenberger parameter |
196 |
|
✗ |
double gamma_L = FelisceParam::instance(this->instanceIndex()).NS_gamma; |
197 |
|
|
|
198 |
|
|
// time derivative coefficient |
199 |
|
✗ |
double rhof = FelisceParam::instance(this->instanceIndex()).density; |
200 |
|
✗ |
double coef_time = rhof/m_fstransient->timeStep; |
201 |
|
|
|
202 |
|
|
// defining the block matrix of the system |
203 |
|
✗ |
if (flagMatrixRHS == FlagMatrixRHS::matrix_and_rhs || flagMatrixRHS == FlagMatrixRHS::only_matrix) { |
204 |
|
|
/// div stabilizing term acting on u, int \nabla \cdot u \nabla \cdot v |
205 |
|
✗ |
m_elementMat[0]->div_phi_j_div_phi_i(gamma_div, *m_feVel, 0, 0, this->dimension()); // block (0,0) |
206 |
|
|
// data term, int_\Omega uv. to this block we will also add the cip stabilization |
207 |
|
✗ |
m_elementMat[0]->phi_i_phi_j(gamma_M, *m_feVel, 0, 0, this->dimension()); // block (0,0) |
208 |
|
|
|
209 |
|
|
// stabilizing term acting on p |
210 |
|
✗ |
m_elementMat[0]->grad_phi_i_grad_phi_j(gamma_p, *m_fePres, this->dimension(), this->dimension(), 1); // block (1,1) |
211 |
|
|
|
212 |
|
|
// time discretization for z |
213 |
|
✗ |
m_elementMat[0]->phi_i_phi_j(coef_time, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
214 |
|
|
|
215 |
|
|
// Stokes blocks for z and y, int_\Omega \epsilon v : \epsilon z or int_\Omega \nabla v : \nabla z |
216 |
|
✗ |
if (m_useSymmetricStress) |
217 |
|
✗ |
m_elementMat[0]->eps_phi_i_eps_phi_j(2*m_viscosity, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
218 |
|
|
else |
219 |
|
✗ |
m_elementMat[0]->grad_phi_i_grad_phi_j(m_viscosity, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
220 |
|
|
// \int y \nabla \cdot v and -\int q \nabla \cdot z |
221 |
|
✗ |
m_elementMat[0]->psi_j_div_phi_i(1., *m_feVel, *m_fePres, 0, 1+2*this->dimension()); // block (0,3) |
222 |
|
✗ |
m_elementMat[0]->psi_i_div_phi_j(-1., *m_feVel, *m_fePres, this->dimension(), 1+this->dimension()); // block (1,2) |
223 |
|
|
|
224 |
|
|
// time discretization for u |
225 |
|
✗ |
m_elementMat[0]->phi_i_phi_j(coef_time, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
226 |
|
|
|
227 |
|
|
// Stokes blocks for u and p, int_\Omega \epsilon u : \epsilon v or int_\Omega \nabla u : \nabla v |
228 |
|
✗ |
if (m_useSymmetricStress) |
229 |
|
✗ |
m_elementMat[0]->eps_phi_i_eps_phi_j(2*m_viscosity, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
230 |
|
|
else |
231 |
|
✗ |
m_elementMat[0]->grad_phi_i_grad_phi_j(m_viscosity, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
232 |
|
|
// -int p \nabla \cdot w and \int x \nabla \cdot u |
233 |
|
✗ |
m_elementMat[0]->psi_j_div_phi_i(-1., *m_feVel, *m_fePres, 1+this->dimension(), this->dimension()); // block (2,1) |
234 |
|
✗ |
m_elementMat[0]->psi_i_div_phi_j(1., *m_feVel, *m_fePres, 1+2*this->dimension(), 0); // block (3,0) |
235 |
|
|
|
236 |
|
|
// stabilizing term acting on z, -int_\Omega grad z : grad w |
237 |
|
✗ |
m_elementMat[0]->grad_phi_i_grad_phi_j(-gamma_u_star, *m_feVel, 1+this->dimension(), 1+this->dimension(), this->dimension()); // block (2,2) |
238 |
|
|
// stabilizing term acting on y, -int_\Omega y * x |
239 |
|
✗ |
m_elementMat[0]->phi_i_phi_j(-gamma_p_star, *m_fePres, 1+2*this->dimension(), 1+2*this->dimension(), 1); // block (3,3) |
240 |
|
|
|
241 |
|
|
// ALE |
242 |
|
✗ |
if(FelisceParam::instance(this->instanceIndex()).useALEformulation) { |
243 |
|
|
// m_elemFielVelMesh -> -w^{n} |
244 |
|
✗ |
m_elemFieldVelMesh.setValue(externalVec(0), *m_feVel, iel, m_iDisplacement, m_externalAO[0], *m_externalDof[0]); |
245 |
|
|
// ALE term for u |
246 |
|
✗ |
m_elementMat[0]->div_u_phi_j_phi_i(rhof, m_elemFieldVelMesh, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
247 |
|
|
// ALE term for z |
248 |
|
✗ |
m_elementMat[0]->div_u_phi_j_phi_i(rhof, m_elemFieldVelMesh, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
249 |
|
|
|
250 |
|
|
|
251 |
|
✗ |
m_beta.setValue(previousVelocity(), *m_feVel, iel, m_iVelocity, *m_aoPreviousVelocity, *m_dofPreviousVelocity); |
252 |
|
|
// Construction of m_beta = u^{n-1} - w^{n}; Recall that externalVec(0) = - w^{n} |
253 |
|
✗ |
m_beta.get_val() = m_beta.get_val() + m_elemFieldVelMesh.get_val(); |
254 |
|
|
|
255 |
|
✗ |
if(FelisceParam::instance(this->instanceIndex()).NSequationFlag) { // for Navier-Stokes |
256 |
|
|
// Temam's stabilization trick |
257 |
|
✗ |
m_elemFieldAdv.setValue(previousVelocity(), *m_feVel, iel, m_iVelocity, aoPreviousVelocity(), dofPreviousVelocity()); |
258 |
|
✗ |
m_elementMat[0]->div_u_phi_j_phi_i(0.5*rhof, m_elemFieldAdv, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
259 |
|
✗ |
m_elementMat[0]->div_u_phi_j_phi_i(0.5*rhof, m_elemFieldAdv, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
260 |
|
|
|
261 |
|
|
|
262 |
|
|
//convective term for u |
263 |
|
✗ |
m_elementMat[0]->u_grad_phi_j_phi_i(rhof, m_beta, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
264 |
|
|
//convective term for z |
265 |
|
✗ |
m_elementMat[0]->u_grad_phi_j_phi_i(rhof, m_beta, *m_feVel, 0, 1+this->dimension(), this->dimension(), true); // block (0,2) |
266 |
|
|
} |
267 |
|
|
else { // for Stokes |
268 |
|
|
//convective ALE term for u |
269 |
|
✗ |
m_elementMat[0]->u_grad_phi_j_phi_i(rhof, m_elemFieldVelMesh, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
270 |
|
|
//convective ALE term for z |
271 |
|
✗ |
m_elementMat[0]->u_grad_phi_j_phi_i(rhof, m_elemFieldVelMesh, *m_feVel, 0, 1+this->dimension(), this->dimension(), true); // block (0,2) |
272 |
|
|
} |
273 |
|
|
} |
274 |
|
|
else { |
275 |
|
✗ |
if(FelisceParam::instance(this->instanceIndex()).NSequationFlag) { // for Navier-Stokes, no ALE |
276 |
|
|
// Temam's stabilization trick |
277 |
|
✗ |
m_elemFieldAdv.setValue(previousVelocity(), *m_feVel, iel, m_iVelocity, aoPreviousVelocity(), dofPreviousVelocity()); |
278 |
|
|
|
279 |
|
✗ |
m_elementMat[0]->div_u_phi_j_phi_i(0.5*rhof, m_elemFieldAdv, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
280 |
|
✗ |
m_elementMat[0]->div_u_phi_j_phi_i(0.5*rhof, m_elemFieldAdv, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
281 |
|
|
|
282 |
|
|
// linearized convective term |
283 |
|
✗ |
m_elementMat[0]->u_grad_phi_j_phi_i(rhof, m_elemFieldAdv, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
284 |
|
✗ |
m_elementMat[0]->u_grad_phi_j_phi_i(rhof, m_elemFieldAdv, *m_feVel, 0,1+this->dimension(),this->dimension(),true); // block (0,2) |
285 |
|
|
} |
286 |
|
|
} |
287 |
|
|
|
288 |
|
|
// Luenberger terms |
289 |
|
|
// int_\Omega v * z |
290 |
|
✗ |
m_elementMat[0]->phi_i_phi_j(gamma_L, *m_feVel, 0, 1+this->dimension(), this->dimension()); // block (0,2) |
291 |
|
|
// int_\Omega u * w |
292 |
|
✗ |
m_elementMat[0]->phi_i_phi_j(gamma_L, *m_feVel, 1+this->dimension(), 0, this->dimension()); // block (2,0) |
293 |
|
|
} |
294 |
|
|
|
295 |
|
|
// rhs terms in the current configuration |
296 |
|
✗ |
if (flagMatrixRHS == FlagMatrixRHS::matrix_and_rhs) { |
297 |
|
|
//rhs data terms, using the measurements |
298 |
|
✗ |
m_dataTerm.setValue(m_velocityData, *m_feVel, iel, m_iVelocity, m_ao, dof()); |
299 |
|
|
|
300 |
|
|
// gamma_M int_\Omega u_M * v |
301 |
|
✗ |
m_elementVector[0]->source(gamma_M, *m_feVel, m_dataTerm, 0, this->dimension()); |
302 |
|
|
|
303 |
|
|
// gamma_L int_\Omega u_M * w |
304 |
|
✗ |
m_elementVector[0]->source(gamma_L, *m_feVel, m_dataTerm, 1+this->dimension(), this->dimension()); // Luenberger term |
305 |
|
|
} |
306 |
|
|
|
307 |
|
|
// rhs terms in the previous configuration |
308 |
|
✗ |
if (flagMatrixRHS == FlagMatrixRHS::only_rhs) { |
309 |
|
|
// rhof/tau int_\Omega u^{n-1} * w |
310 |
|
✗ |
m_sourceTerm.setValue(previousVelocity(), *m_feVel, iel, m_iVelocity, aoPreviousVelocity(), dofPreviousVelocity()); |
311 |
|
|
|
312 |
|
✗ |
m_elementVector[0]->source(coef_time, *m_feVel, m_sourceTerm, 1+this->dimension(), this->dimension()); |
313 |
|
|
// rhs f term, zero |
314 |
|
|
} |
315 |
|
|
} |
316 |
|
|
|
317 |
|
✗ |
void LinearProblemNSContinuation::userChangePattern(int numProc, int rankProc) { |
318 |
|
|
IGNORE_UNUSED_ARGUMENT(numProc); |
319 |
|
|
IGNORE_UNUSED_ARGUMENT(rankProc); |
320 |
|
|
// compute initial (uniform) repartition of dof on processes |
321 |
|
✗ |
felInt numDofByProc = m_numDof/MpiInfo::numProc(); |
322 |
|
✗ |
std::vector<felInt> numDofPerProc(MpiInfo::numProc()); |
323 |
|
✗ |
for(felInt i=0; i<MpiInfo::numProc(); ++i) { |
324 |
|
✗ |
if(i == MpiInfo::numProc() - 1) |
325 |
|
✗ |
numDofPerProc[i] = m_numDof - i*numDofByProc; |
326 |
|
|
else |
327 |
|
✗ |
numDofPerProc[i] = numDofByProc; |
328 |
|
|
} |
329 |
|
|
|
330 |
|
✗ |
felInt shiftDof = 0; |
331 |
|
✗ |
for(felInt i=0; i<MpiInfo::rankProc(); ++i) |
332 |
|
✗ |
shiftDof += numDofPerProc[i]; |
333 |
|
|
|
334 |
|
✗ |
std::vector<felInt> rankDof(m_numDof, -1); |
335 |
|
✗ |
for(felInt i=0; i<numDofPerProc[MpiInfo::rankProc()]; ++i) |
336 |
|
✗ |
rankDof[i + shiftDof] = MpiInfo::rankProc(); |
337 |
|
|
|
338 |
|
|
// build the edges or faces depending on the dimension (for the global mesh) |
339 |
|
|
// In this function, "faces" refers to either edges or faces. |
340 |
|
✗ |
if(dimension() == 2) |
341 |
|
✗ |
m_mesh[m_currentMesh]->buildEdges(); |
342 |
|
|
else |
343 |
|
✗ |
m_mesh[m_currentMesh]->buildFaces(); |
344 |
|
|
|
345 |
|
|
// variables |
346 |
|
|
GeometricMeshRegion::ElementType eltType, eltTypeOpp; |
347 |
|
|
felInt idVar1, idVar2; |
348 |
|
|
felInt node1, node2; |
349 |
|
|
felInt numFacesPerElement, numEltPerLabel; |
350 |
|
|
felInt ielSupport, jelSupport; |
351 |
|
|
felInt ielCurrentGeo, ielOppositeGeo; |
352 |
|
✗ |
std::vector<felInt> vecSupport, vecSupportOpposite; |
353 |
|
✗ |
std::vector<felInt> numElement(m_mesh[m_currentMesh]->m_numTypesOfElement, 0); |
354 |
|
✗ |
std::vector< std::set<felInt> > nodesNeighborhood(m_numDof); |
355 |
|
|
|
356 |
|
|
|
357 |
|
|
// zeroth loop on the unknown of the linear problem |
358 |
|
✗ |
for (size_t iUnknown1 = 0; iUnknown1 < m_listUnknown.size(); ++iUnknown1) { |
359 |
|
✗ |
idVar1 = m_listUnknown.idVariable(iUnknown1); |
360 |
|
|
|
361 |
|
✗ |
ielCurrentGeo = 0; |
362 |
|
|
|
363 |
|
|
// first loop on element type |
364 |
|
✗ |
for (size_t i=0; i<m_mesh[m_currentMesh]->bagElementTypeDomain().size(); ++i) { |
365 |
|
✗ |
eltType = m_mesh[m_currentMesh]->bagElementTypeDomain()[i]; |
366 |
|
✗ |
const GeoElement* geoEle = GeometricMeshRegion::eltEnumToFelNameGeoEle[eltType].second; |
367 |
|
✗ |
numElement[eltType] = 0; |
368 |
|
✗ |
numFacesPerElement = geoEle->numBdEle(); |
369 |
|
|
|
370 |
|
|
// second loop on region of the mesh. |
371 |
|
✗ |
for(GeometricMeshRegion::IntRefToBegEndIndex_type::const_iterator itRef =m_mesh[m_currentMesh]->intRefToBegEndMaps[eltType].begin(); itRef !=m_mesh[m_currentMesh]->intRefToBegEndMaps[eltType].end(); itRef++) { |
372 |
|
✗ |
numEltPerLabel = itRef->second.second; |
373 |
|
|
|
374 |
|
|
// Third loop on element |
375 |
|
✗ |
for (felInt iel = 0; iel < numEltPerLabel; iel++) { |
376 |
|
✗ |
m_supportDofUnknown[iUnknown1].getIdElementSupport(eltType, numElement[eltType], vecSupport); |
377 |
|
|
|
378 |
|
✗ |
for(size_t ielSup=0; ielSup<vecSupport.size(); ++ielSup) { |
379 |
|
✗ |
ielSupport = vecSupport[ielSup]; |
380 |
|
|
|
381 |
|
|
// for all support dof in this support element |
382 |
|
✗ |
for (felInt iSup = 0; iSup < m_supportDofUnknown[iUnknown1].getNumSupportDof(ielSupport); ++iSup) { |
383 |
|
|
// for all component of the unknown |
384 |
|
✗ |
for (size_t iComp = 0; iComp < m_listVariable[idVar1].numComponent(); iComp++) { |
385 |
|
|
// get the global id of the support dof |
386 |
|
✗ |
dof().loc2glob(ielSupport, iSup, idVar1, iComp, node1); |
387 |
|
|
|
388 |
|
|
// if this node is on this process |
389 |
|
✗ |
if(rankDof[node1] == MpiInfo::rankProc()) { |
390 |
|
|
// loop over the boundaries of the element |
391 |
|
✗ |
for(felInt iface=0; iface < numFacesPerElement; ++iface) { |
392 |
|
|
// check if this face is an inner face or a boundary |
393 |
|
✗ |
ielOppositeGeo = ielCurrentGeo; |
394 |
|
✗ |
eltTypeOpp = eltType; |
395 |
|
✗ |
bool isAnInnerFace =m_mesh[m_currentMesh]->getAdjElement(eltTypeOpp, ielOppositeGeo, iface); |
396 |
|
|
|
397 |
|
✗ |
if(isAnInnerFace) { |
398 |
|
|
// for all unknown |
399 |
|
✗ |
for (size_t iUnknown2 = 0; iUnknown2 < m_listUnknown.size(); ++iUnknown2) { |
400 |
|
✗ |
idVar2 = m_listUnknown.idVariable(iUnknown2); |
401 |
|
|
|
402 |
|
|
// get the support element of the opposite element |
403 |
|
✗ |
m_supportDofUnknown[iUnknown2].getIdElementSupport(ielOppositeGeo, vecSupportOpposite); |
404 |
|
|
|
405 |
|
|
// for all support element |
406 |
|
✗ |
for(size_t jelSup=0; jelSup<vecSupportOpposite.size(); ++jelSup) { |
407 |
|
✗ |
jelSupport = vecSupportOpposite[jelSup]; |
408 |
|
|
|
409 |
|
|
// for all support dof in this support element |
410 |
|
✗ |
for (felInt jSup = 0; jSup < m_supportDofUnknown[iUnknown2].getNumSupportDof(jelSupport); ++jSup) { |
411 |
|
|
|
412 |
|
|
// for all component of the second unknown |
413 |
|
✗ |
for (size_t jComp = 0; jComp < m_listVariable[idVar2].numComponent(); jComp++) { |
414 |
|
|
|
415 |
|
|
// If the two current components are connected |
416 |
|
✗ |
felInt iConnect = dof().getNumGlobComp(iUnknown1, iComp); |
417 |
|
✗ |
felInt jConnect = dof().getNumGlobComp(iUnknown2, jComp); |
418 |
|
✗ |
if (m_listUnknown.mask()(iConnect, jConnect) > 0) { |
419 |
|
|
// get the global id of the second support dof |
420 |
|
✗ |
dof().loc2glob(jelSupport, jSup, idVar2, jComp, node2); |
421 |
|
|
|
422 |
|
|
// remove diagonal term to define pattern and use it in Parmetis. |
423 |
|
|
// if the two global ids are different, they are neighboors |
424 |
|
✗ |
if(node1 != node2) |
425 |
|
✗ |
nodesNeighborhood[node1].insert(node2); |
426 |
|
|
} |
427 |
|
|
} |
428 |
|
|
} |
429 |
|
|
} |
430 |
|
|
} |
431 |
|
|
} |
432 |
|
|
} |
433 |
|
|
} |
434 |
|
|
} |
435 |
|
|
} |
436 |
|
|
} |
437 |
|
✗ |
++ielCurrentGeo; |
438 |
|
✗ |
++numElement[eltType]; |
439 |
|
|
} |
440 |
|
|
} |
441 |
|
|
} |
442 |
|
|
} |
443 |
|
|
|
444 |
|
|
// Store the pattern in CSR style |
445 |
|
✗ |
std::vector<felInt> iCSR, jCSR; |
446 |
|
✗ |
felInt dofSize = 0; |
447 |
|
✗ |
felInt cptDof = 0; |
448 |
|
|
felInt pos; |
449 |
|
|
|
450 |
|
✗ |
iCSR.resize(dof().pattern().numRows() + 1, 0); |
451 |
|
✗ |
for(size_t iNode=0; iNode<nodesNeighborhood.size(); ++iNode) |
452 |
|
✗ |
dofSize += nodesNeighborhood[iNode].size(); |
453 |
|
|
|
454 |
|
✗ |
jCSR.resize(dofSize, 0); |
455 |
|
✗ |
for(size_t iNode=0; iNode<nodesNeighborhood.size(); ++iNode) { |
456 |
|
✗ |
if(rankDof[iNode] == MpiInfo::rankProc()) { |
457 |
|
✗ |
iCSR[cptDof + 1] = iCSR[cptDof] + nodesNeighborhood[iNode].size(); |
458 |
|
✗ |
pos = 0; |
459 |
|
✗ |
for(std::set<felInt>::iterator it=nodesNeighborhood[iNode].begin(); it != nodesNeighborhood[iNode].end(); ++it) { |
460 |
|
✗ |
jCSR[iCSR[cptDof] + pos] = *it; |
461 |
|
✗ |
++pos; |
462 |
|
|
} |
463 |
|
✗ |
++cptDof; |
464 |
|
|
} |
465 |
|
|
} |
466 |
|
|
|
467 |
|
|
// Now, call merge pattern |
468 |
|
✗ |
dof().mergeGlobalPattern(iCSR, jCSR); |
469 |
|
|
//std::cout<<"After global merge"<<std::endl; |
470 |
|
|
//dof().pattern().print(1, std::cout); |
471 |
|
|
} |
472 |
|
|
|
473 |
|
|
// compared to the function linearProblemNS::assembleFaceOrientedStabilization, here we use different blocks and only the jumps in the velocity |
474 |
|
✗ |
void LinearProblemNSContinuation::assembleCIPStabilization() { |
475 |
|
|
// definition of variables |
476 |
|
✗ |
std::pair<bool, GeometricMeshRegion::ElementType> adjElt; |
477 |
|
|
GeometricMeshRegion::ElementType eltType; // Type of element |
478 |
|
|
|
479 |
|
✗ |
felInt ielCurrentLocalGeo = 0; // local geometric id of the current element |
480 |
|
✗ |
felInt ielCurrentGlobalGeo = 0; // global geometric id of the current element |
481 |
|
✗ |
felInt ielOppositeGlobalGeo = 0; // global geometric id of the opposite element |
482 |
|
|
felInt numFacesPerType; // number of faces of an element of a given type |
483 |
|
|
felInt numPointPerElt; // number of vertices by element |
484 |
|
|
felInt ielCurrentGlobalSup; // global support element id of the current element |
485 |
|
|
felInt ielOppositeGlobalSup; // global support element id of the opposite element |
486 |
|
|
|
487 |
|
✗ |
std::vector<felInt> idOfFacesCurrent; // ids of the current element edges |
488 |
|
✗ |
std::vector<felInt> idOfFacesOpposite; // ids of the opposite element edges |
489 |
|
✗ |
std::vector<felInt> currentElemIdPoint; // ids of the vertices of the current element |
490 |
|
✗ |
std::vector<felInt> oppositeElemIdPoint; // ids of the vertices of the opposite element |
491 |
|
✗ |
std::vector<Point*> currentElemPoint; // point coordinates of the current element vertices |
492 |
|
✗ |
std::vector<Point*> oppositeElemPoint; // point coordinates of the opposite element vertices |
493 |
|
|
|
494 |
|
✗ |
FlagMatrixRHS flag = FlagMatrixRHS::only_matrix; // flag to only assemble the matrix |
495 |
|
|
felInt idFaceToDo; |
496 |
|
✗ |
bool allDone = false; |
497 |
|
|
|
498 |
|
✗ |
ElementField elemFieldAdvFace; |
499 |
|
|
|
500 |
|
|
felInt rank; |
501 |
|
✗ |
MPI_Comm_rank(PETSC_COMM_WORLD, &rank); |
502 |
|
|
|
503 |
|
|
// bag element type vector |
504 |
|
✗ |
const std::vector<ElementType>& bagElementTypeDomain = m_meshLocal[m_currentMesh]->bagElementTypeDomain(); |
505 |
|
✗ |
eltType = bagElementTypeDomain[0]; |
506 |
|
|
|
507 |
|
|
|
508 |
|
|
// Finite Element With Bd for the opposite element |
509 |
|
|
CurrentFiniteElementWithBd* ptmp; |
510 |
|
|
CurrentFiniteElementWithBd* oppositeFEWithBd; |
511 |
|
|
|
512 |
|
✗ |
m_velocity = &m_listVariable[m_listVariable.getVariableIdList(velocity)]; |
513 |
|
|
|
514 |
|
✗ |
const GeoElement *geoEle = m_mesh[m_currentMesh]->eltEnumToFelNameGeoEle[eltType].second; |
515 |
|
✗ |
felInt typeOfFEVel = m_velocity->finiteElementType(); |
516 |
|
|
|
517 |
|
✗ |
const RefElement *refEleVel = geoEle->defineFiniteEle(eltType, typeOfFEVel, *m_mesh[m_currentMesh]); |
518 |
|
✗ |
oppositeFEWithBd = new CurrentFiniteElementWithBd(*refEleVel, *geoEle, m_velocity->degreeOfExactness(), m_velocity->degreeOfExactness()); |
519 |
|
|
|
520 |
|
|
// initializing variables |
521 |
|
✗ |
numPointPerElt = m_meshLocal[m_currentMesh]->m_numPointsPerElt[eltType]; |
522 |
|
✗ |
numFacesPerType = GeometricMeshRegion::eltEnumToFelNameGeoEle[eltType].second->numBdEle(); |
523 |
|
|
|
524 |
|
|
// resize of vectors |
525 |
|
✗ |
idOfFacesCurrent.resize(numFacesPerType, -1); |
526 |
|
✗ |
idOfFacesOpposite.resize(numFacesPerType, -1); |
527 |
|
✗ |
currentElemIdPoint.resize(numPointPerElt, -1); |
528 |
|
✗ |
oppositeElemIdPoint.resize(numPointPerElt, -1); |
529 |
|
✗ |
currentElemPoint.resize(numPointPerElt, NULL); |
530 |
|
✗ |
oppositeElemPoint.resize(numPointPerElt, NULL); |
531 |
|
|
|
532 |
|
|
// define finite element |
533 |
|
✗ |
defineFiniteElement(eltType); |
534 |
|
✗ |
initElementArray(); |
535 |
|
✗ |
defineCurrentFiniteElementWithBd(eltType); |
536 |
|
|
|
537 |
|
|
// allocate arrays for assembling the matrix |
538 |
|
✗ |
allocateArrayForAssembleMatrixRHS(flag); |
539 |
|
|
|
540 |
|
|
// init variables |
541 |
|
✗ |
initPerElementType(eltType, flag); |
542 |
|
✗ |
CurrentFiniteElementWithBd* feWithBd = m_listCurrentFiniteElementWithBd[m_iVelocity]; |
543 |
|
|
|
544 |
|
✗ |
CurrentFiniteElementWithBd* firstCurrentFe = feWithBd; |
545 |
|
✗ |
CurrentFiniteElementWithBd* firstOppositeFe = oppositeFEWithBd; |
546 |
|
|
|
547 |
|
|
// get informations on the current element |
548 |
|
✗ |
setElemPoint(eltType, 0, currentElemPoint, currentElemIdPoint, &ielCurrentGlobalSup); |
549 |
|
|
|
550 |
|
|
// update the finite elements |
551 |
|
✗ |
feWithBd->updateFirstDeriv(0, currentElemPoint); |
552 |
|
✗ |
feWithBd->updateBdMeasNormal(); |
553 |
|
|
|
554 |
|
|
// get the global id of the first geometric element |
555 |
|
✗ |
ISLocalToGlobalMappingApply(m_mappingElem[m_currentMesh], 1, &ielCurrentLocalGeo, &ielCurrentGlobalGeo); |
556 |
|
|
|
557 |
|
|
// the map to remember what contribution have already been computed |
558 |
|
✗ |
std::map<felInt, std::vector<bool> > contribDone; |
559 |
|
✗ |
addNewFaceOrientedContributor(numFacesPerType, ielCurrentGlobalGeo, contribDone[ielCurrentGlobalGeo]); |
560 |
|
|
|
561 |
|
|
// build the stack to know what is the next element |
562 |
|
✗ |
std::stack<felInt> nextInStack; |
563 |
|
✗ |
nextInStack.push(ielCurrentGlobalGeo); |
564 |
|
|
|
565 |
|
|
// get all the faces of the element |
566 |
|
✗ |
if(dimension() == 2) |
567 |
|
✗ |
m_mesh[m_currentMesh]->getAllEdgeOfElement(ielCurrentGlobalGeo, idOfFacesCurrent); |
568 |
|
|
else |
569 |
|
✗ |
m_mesh[m_currentMesh]->getAllFaceOfElement(ielCurrentGlobalGeo, idOfFacesCurrent); |
570 |
|
|
|
571 |
|
|
// loop over all the element (snake style) |
572 |
|
✗ |
while(!nextInStack.empty()) { |
573 |
|
|
// check the faces and use the first one that is an inner face and that has not been done yet. |
574 |
|
✗ |
idFaceToDo = -1; |
575 |
|
✗ |
allDone = true; |
576 |
|
✗ |
for(size_t iface=0; iface<contribDone[ielCurrentGlobalGeo].size(); ++iface) { |
577 |
|
✗ |
if(!contribDone[ielCurrentGlobalGeo][iface]) { |
578 |
|
|
// This is an inner face, check if the contribution has already been computed or not |
579 |
|
✗ |
if(idFaceToDo == -1) { |
580 |
|
✗ |
idFaceToDo = iface; |
581 |
|
|
} else { |
582 |
|
✗ |
allDone = false; |
583 |
|
|
} |
584 |
|
|
} |
585 |
|
|
} |
586 |
|
|
|
587 |
|
|
// update the stack |
588 |
|
✗ |
if(!allDone) |
589 |
|
✗ |
nextInStack.push(ielCurrentGlobalGeo); |
590 |
|
|
|
591 |
|
✗ |
if(nextInStack.top() == ielCurrentGlobalGeo && allDone) |
592 |
|
✗ |
nextInStack.pop(); |
593 |
|
|
|
594 |
|
|
|
595 |
|
|
// assemble terms |
596 |
|
✗ |
if(idFaceToDo != -1) { |
597 |
|
|
// get the opposite id of the element |
598 |
|
✗ |
ielOppositeGlobalGeo = ielCurrentGlobalGeo; |
599 |
|
✗ |
adjElt =m_mesh[m_currentMesh]->getAdjElement(ielOppositeGlobalGeo, idFaceToDo); |
600 |
|
|
|
601 |
|
|
// get the type of the opposite element |
602 |
|
|
// eltTypeOpp = adjElt.second; |
603 |
|
|
|
604 |
|
|
// update the opposite finite element and set ielOppositeGlobalSup |
605 |
|
✗ |
updateFaceOrientedFEWithBd(oppositeFEWithBd, idOfFacesOpposite, numPointPerElt, ielOppositeGlobalGeo, ielOppositeGlobalSup); |
606 |
|
|
|
607 |
|
|
// find the local id of the edge in the opposite element |
608 |
|
✗ |
felInt localIdFaceOpposite = -1; |
609 |
|
✗ |
for(size_t jface=0; jface<idOfFacesOpposite.size(); ++jface) { |
610 |
|
✗ |
if(idOfFacesCurrent[idFaceToDo] == idOfFacesOpposite[jface]) { |
611 |
|
✗ |
localIdFaceOpposite = jface; |
612 |
|
|
} |
613 |
|
|
} |
614 |
|
|
|
615 |
|
|
// compute coefficients |
616 |
|
✗ |
CurvilinearFiniteElement* curvFe = &feWithBd->bdEle(idFaceToDo); |
617 |
|
|
//CurvilinearFiniteElement* curvFeOP = &oppositeFEWithBd->bdEle(localIdFaceOpposite); |
618 |
|
✗ |
double hK = curvFe->diameter(); //0.5 * (feWithBd->diameter() + oppositeFEWithBd->diameter()); |
619 |
|
|
//cout << hK - curvFeOP->diameter() << endl; |
620 |
|
✗ |
double gamma_u = FelisceParam::instance(this->instanceIndex()).referenceValue1; // cip coefficient |
621 |
|
✗ |
double gamma = hK*gamma_u; |
622 |
|
|
|
623 |
|
|
// We can now compute all the integrals. |
624 |
|
|
// current element for phi_i and phi_j |
625 |
|
✗ |
m_elementMat[0]->zero(); |
626 |
|
✗ |
m_elementMat[0]->grad_phi_j_dot_n_grad_psi_i_dot_n(gamma, *feWithBd, *feWithBd, idFaceToDo, idFaceToDo, 0, 0, this->dimension()); // block (0,0) |
627 |
|
✗ |
setValueMatrixRHS(ielCurrentGlobalSup, ielCurrentGlobalSup, flag); |
628 |
|
|
|
629 |
|
|
// current element for phi_i and opposite element for phi_j |
630 |
|
✗ |
m_elementMat[0]->zero(); |
631 |
|
✗ |
m_elementMat[0]->grad_phi_j_dot_n_grad_psi_i_dot_n(gamma, *oppositeFEWithBd, *feWithBd, localIdFaceOpposite, idFaceToDo, 0, 0, this->dimension()); // block (0,0) |
632 |
|
✗ |
setValueMatrixRHS(ielCurrentGlobalSup, ielOppositeGlobalSup, flag); |
633 |
|
|
|
634 |
|
|
// opposite element for phi_i and current element for phi_j |
635 |
|
✗ |
m_elementMat[0]->zero(); |
636 |
|
✗ |
m_elementMat[0]->grad_phi_j_dot_n_grad_psi_i_dot_n(gamma, *feWithBd, *oppositeFEWithBd, idFaceToDo, localIdFaceOpposite, 0, 0, this->dimension()); // block (0,0) |
637 |
|
✗ |
setValueMatrixRHS(ielOppositeGlobalSup, ielCurrentGlobalSup, flag); |
638 |
|
|
|
639 |
|
|
// opposite element for phi_i and phi_j |
640 |
|
✗ |
m_elementMat[0]->zero(); |
641 |
|
✗ |
m_elementMat[0]->grad_phi_j_dot_n_grad_psi_i_dot_n(gamma, *oppositeFEWithBd, *oppositeFEWithBd, localIdFaceOpposite, localIdFaceOpposite, 0, 0, this->dimension()); // block (0,0) |
642 |
|
✗ |
setValueMatrixRHS(ielOppositeGlobalSup, ielOppositeGlobalSup, flag); |
643 |
|
|
|
644 |
|
|
// update the map to say that this face is done |
645 |
|
✗ |
contribDone[ielCurrentGlobalGeo][idFaceToDo] = true; |
646 |
|
|
|
647 |
|
|
// check if the opposite element is on this proc |
648 |
|
✗ |
if(m_eltPart[m_currentMesh][ielOppositeGlobalGeo] == rank) { |
649 |
|
✗ |
if(contribDone.find(ielOppositeGlobalGeo) == contribDone.end()) { |
650 |
|
|
// not found in the map, add it and initialize it. |
651 |
|
✗ |
addNewFaceOrientedContributor(numFacesPerType, ielOppositeGlobalGeo, contribDone[ielOppositeGlobalGeo]); |
652 |
|
|
} |
653 |
|
✗ |
contribDone[ielOppositeGlobalGeo][localIdFaceOpposite] = true; |
654 |
|
|
|
655 |
|
|
// update the finite element as the opposite finite element. |
656 |
|
✗ |
ptmp = feWithBd; |
657 |
|
✗ |
feWithBd = oppositeFEWithBd; |
658 |
|
✗ |
oppositeFEWithBd = ptmp; |
659 |
|
|
|
660 |
|
✗ |
ielCurrentGlobalGeo = ielOppositeGlobalGeo; |
661 |
|
✗ |
ielCurrentGlobalSup = ielOppositeGlobalSup; |
662 |
|
✗ |
idOfFacesCurrent = idOfFacesOpposite; |
663 |
|
|
} |
664 |
|
|
else { |
665 |
|
|
// not on this rank, take the next one in the stack |
666 |
|
✗ |
if(!nextInStack.empty()) { |
667 |
|
✗ |
ielCurrentGlobalGeo = nextInStack.top(); |
668 |
|
✗ |
updateFaceOrientedFEWithBd(feWithBd, idOfFacesCurrent, numPointPerElt, ielCurrentGlobalGeo, ielCurrentGlobalSup); |
669 |
|
|
} |
670 |
|
|
} |
671 |
|
|
} |
672 |
|
|
else { |
673 |
|
|
// All contribution have been already computed on this element |
674 |
|
|
// take the next element in the stack |
675 |
|
✗ |
if(!nextInStack.empty()) { |
676 |
|
✗ |
ielCurrentGlobalGeo = nextInStack.top(); |
677 |
|
✗ |
updateFaceOrientedFEWithBd(feWithBd, idOfFacesCurrent, numPointPerElt, ielCurrentGlobalGeo, ielCurrentGlobalSup); |
678 |
|
|
} |
679 |
|
|
} |
680 |
|
|
} |
681 |
|
|
|
682 |
|
|
// desallocate array use for assemble with petsc. |
683 |
|
✗ |
desallocateArrayForAssembleMatrixRHS(flag); |
684 |
|
|
|
685 |
|
|
// desallocate opposite finite elements |
686 |
|
✗ |
feWithBd = firstCurrentFe; |
687 |
|
✗ |
delete firstOppositeFe; |
688 |
|
|
} |
689 |
|
|
|
690 |
|
✗ |
void LinearProblemNSContinuation::addNewFaceOrientedContributor(felInt size, felInt idElt, std::vector<bool>& vec) { |
691 |
|
|
felInt ielTmp; |
692 |
|
✗ |
std::pair<bool, GeometricMeshRegion::ElementType> adjElt; |
693 |
|
|
|
694 |
|
✗ |
vec.resize(size, false); |
695 |
|
✗ |
for(felInt iface=0; iface<size; ++iface) { |
696 |
|
|
// check if this face is an inner face or a boundary |
697 |
|
✗ |
ielTmp = idElt; |
698 |
|
✗ |
adjElt =m_mesh[m_currentMesh]->getAdjElement(ielTmp, iface); |
699 |
|
|
|
700 |
|
✗ |
if(!adjElt.first) { |
701 |
|
|
// not an inner face, set the contribution to done. |
702 |
|
✗ |
vec[iface] = true; |
703 |
|
|
} |
704 |
|
|
} |
705 |
|
|
} |
706 |
|
|
|
707 |
|
✗ |
void LinearProblemNSContinuation::updateFaceOrientedFEWithBd(CurrentFiniteElementWithBd* fe, std::vector<felInt>& idOfFaces, felInt numPoints, felInt idElt, felInt& idSup) { |
708 |
|
✗ |
std::vector<felInt> elemIdPoint(numPoints, -1); |
709 |
|
✗ |
std::vector<Point*> elemPoint(numPoints, NULL); |
710 |
|
|
|
711 |
|
|
// get information on the opposite element |
712 |
|
|
// same as the function "setElemPoint" but here ielOppositeGlobalGeo is global |
713 |
|
|
// we assume that the supportDofMesh is the same for all unknown (like in setElemPoint) |
714 |
|
✗ |
m_supportDofUnknown[0].getIdElementSupport(idElt, idSup); |
715 |
|
✗ |
m_mesh[m_currentMesh]->getOneElement(idElt, elemIdPoint); |
716 |
|
✗ |
for (felInt iPoint=0; iPoint<numPoints; ++iPoint) |
717 |
|
✗ |
elemPoint[iPoint] = &(m_mesh[m_currentMesh]->listPoints()[elemIdPoint[iPoint]]); |
718 |
|
|
|
719 |
|
|
// update the finite elements |
720 |
|
✗ |
fe->updateFirstDeriv(0, elemPoint); |
721 |
|
✗ |
fe->updateBdMeasNormal(); |
722 |
|
|
|
723 |
|
|
// get all the faces of the element |
724 |
|
✗ |
if(dimension() == 2) |
725 |
|
✗ |
m_mesh[m_currentMesh]->getAllEdgeOfElement(idElt, idOfFaces); |
726 |
|
|
else |
727 |
|
✗ |
m_mesh[m_currentMesh]->getAllFaceOfElement(idElt, idOfFaces); |
728 |
|
|
} |
729 |
|
|
} |
730 |
|
|
|