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: |
13 |
|
|
// |
14 |
|
|
|
15 |
|
|
// System includes |
16 |
|
|
|
17 |
|
|
// External includes |
18 |
|
|
|
19 |
|
|
// Project includes |
20 |
|
|
#include "Model/NSFracStepExplicitModel.hpp" |
21 |
|
|
#include "Solver/linearProblemNSFracStepProj.hpp" |
22 |
|
|
|
23 |
|
|
namespace felisce { |
24 |
|
✗ |
NSFracStepExplicitModel::NSFracStepExplicitModel():Model() { |
25 |
|
✗ |
m_lumpedModelBC = nullptr; |
26 |
|
✗ |
m_name = "NSFracStepExplicit"; |
27 |
|
|
} |
28 |
|
|
|
29 |
|
|
|
30 |
|
✗ |
NSFracStepExplicitModel::~NSFracStepExplicitModel() { |
31 |
|
✗ |
if(m_lumpedModelBC) delete m_lumpedModelBC; |
32 |
|
✗ |
m_solExtrapolate.destroy(); |
33 |
|
|
} |
34 |
|
|
|
35 |
|
|
|
36 |
|
✗ |
void NSFracStepExplicitModel::initializeDerivedModel() { |
37 |
|
✗ |
if( FelisceParam::instance().lumpedModelBCLabel.size() ) { |
38 |
|
✗ |
m_lumpedModelBC = new LumpedModelBC(m_fstransient); |
39 |
|
✗ |
static_cast<LinearProblemNSFracStepProj*>(m_linearProblem[2])->initializeLumpedModelBC(m_lumpedModelBC); |
40 |
|
|
} |
41 |
|
|
|
42 |
|
|
//if (FelisceParam::instance().orderBdfFS > 2) |
43 |
|
|
// FEL_ERROR("BDF not yet implemented for order greater than 2 with fractional step model."); |
44 |
|
✗ |
if (FelisceParam::instance().orderBdfFS > 1) |
45 |
|
✗ |
FEL_ERROR("BDF not yet implemented for explicit fractional step model."); |
46 |
|
|
|
47 |
|
✗ |
m_bdfFS.defineOrder(FelisceParam::instance().orderBdfFS); |
48 |
|
✗ |
m_linearProblem[0]->initializeTimeScheme(&m_bdfFS); |
49 |
|
✗ |
m_linearProblem[1]->initializeTimeScheme(&m_bdfFS); |
50 |
|
|
} |
51 |
|
|
|
52 |
|
|
// std::set external vectors for linear problems |
53 |
|
|
// (to be used in assembling routines) |
54 |
|
✗ |
void NSFracStepExplicitModel::setExternalVec() { |
55 |
|
|
// problem 0: |
56 |
|
|
// advection -> rho/tau*(w - un) + un*grad(w) + stab = 0 |
57 |
|
|
// problem 1: |
58 |
|
|
// diffusion -> rho/tau(u,v) + (eps(u),eps(v)) = rho/tau(w,v) - k(grad pn,v) |
59 |
|
|
// problem 2: |
60 |
|
|
// projection -> (grad p,grad q) = -rho/tau (div(u),q) |
61 |
|
✗ |
if (FelisceParam::instance().explicitAdvection) { |
62 |
|
|
// explicit advection needs: u_diffusion, pressure |
63 |
|
✗ |
m_linearProblem[0]->pushBackExternalVec(m_linearProblem[1]->sequentialSolution()); |
64 |
|
✗ |
m_linearProblem[0]->pushBackExternalAO(m_linearProblem[1]->ao()); |
65 |
|
✗ |
m_linearProblem[0]->pushBackExternalDof(m_linearProblem[1]->dof()); |
66 |
|
|
|
67 |
|
|
// diffusion needs: u_advection, pressure (optional) |
68 |
|
✗ |
m_linearProblem[1]->pushBackExternalVec(m_linearProblem[0]->sequentialSolution()); |
69 |
|
✗ |
m_linearProblem[1]->pushBackExternalAO(m_linearProblem[0]->ao()); |
70 |
|
✗ |
m_linearProblem[1]->pushBackExternalDof(m_linearProblem[0]->dof()); |
71 |
|
|
|
72 |
|
✗ |
m_linearProblem[1]->pushBackExternalVec(m_linearProblem[2]->sequentialSolution()); |
73 |
|
✗ |
m_linearProblem[1]->pushBackExternalAO(m_linearProblem[2]->ao()); |
74 |
|
✗ |
m_linearProblem[1]->pushBackExternalDof(m_linearProblem[2]->dof()); |
75 |
|
|
|
76 |
|
|
// projection needs: u_diffusion |
77 |
|
✗ |
m_linearProblem[2]->pushBackExternalVec(m_linearProblem[1]->sequentialSolution()); |
78 |
|
✗ |
m_linearProblem[2]->pushBackExternalAO(m_linearProblem[1]->ao()); |
79 |
|
✗ |
m_linearProblem[2]->pushBackExternalDof(m_linearProblem[1]->dof()); |
80 |
|
|
|
81 |
|
|
} else { |
82 |
|
✗ |
FEL_ERROR(" This model should be used only with explicitAdvection=true"); |
83 |
|
|
} |
84 |
|
|
} |
85 |
|
|
|
86 |
|
|
// --- |
87 |
|
|
// for each problem: assemble matrix and RHS |
88 |
|
|
// first: assemble static matrix (at iteration 0) |
89 |
|
✗ |
void NSFracStepExplicitModel::preAssemblingMatrixRHS(std::size_t iProblem) { |
90 |
|
✗ |
std::cout << " [*] preAssemblingMatrixRHS for problem " << iProblem << " [*] " << std::endl; |
91 |
|
✗ |
m_solExtrapolate.duplicateFrom(m_linearProblem[0]->vector()); |
92 |
|
✗ |
m_solExtrapolate.zeroEntries(); |
93 |
|
|
|
94 |
|
|
// todo : initialize solution for solver. |
95 |
|
|
// m_linearProblem[0]->set_U_0(_U_0); |
96 |
|
✗ |
std::cout << " ... assembleMatrixRHS ... " << std::endl; |
97 |
|
|
// First assembling loop in iteration 0 to build static matrix. |
98 |
|
✗ |
m_linearProblem[iProblem]->assembleMatrixRHS(MpiInfo::rankProc()); |
99 |
|
|
// save static matrix in matrix _A. |
100 |
|
✗ |
std::cout << " ... copyMatrixRHS ... " << std::endl; |
101 |
|
✗ |
m_linearProblem[iProblem]->copyMatrixRHS(); |
102 |
|
|
} |
103 |
|
|
|
104 |
|
|
// first: add dynamic matrix (assembled at each step) |
105 |
|
✗ |
void NSFracStepExplicitModel::postAssemblingMatrixRHS(std::size_t iProblem) { |
106 |
|
✗ |
std::cout << " [*] preAssemblingMatrixRHS for problem " << iProblem << " [*] " << std::endl; |
107 |
|
|
// _Matrix = _A + _Matrix (add static matrix to the dynamic matrix to build |
108 |
|
|
// complete matrix of the system |
109 |
|
✗ |
m_linearProblem[iProblem]->addMatrixRHS(); |
110 |
|
|
} |
111 |
|
|
|
112 |
|
|
// advancing |
113 |
|
✗ |
void NSFracStepExplicitModel::forward() { |
114 |
|
|
// Write solution for postprocessing (if required) |
115 |
|
✗ |
writeSolution(); |
116 |
|
|
// Advance time step. |
117 |
|
✗ |
updateTime(); |
118 |
|
|
// Print time information |
119 |
|
✗ |
printNewTimeIterationBanner(); |
120 |
|
|
|
121 |
|
✗ |
for(std::size_t iLinPb=0; iLinPb<m_linearProblem.size(); iLinPb++) { |
122 |
|
✗ |
if (FelisceParam::verbose()) |
123 |
|
✗ |
PetscPrintf(MpiInfo::petscComm(), "NSFracStepExplicit: linear problem #%ld \n",static_cast<unsigned long>(iLinPb)); |
124 |
|
|
|
125 |
|
|
// initialize BDF for advection |
126 |
|
✗ |
if ( iLinPb==0 ) { |
127 |
|
✗ |
if ( m_fstransient->iteration == 1 ) { |
128 |
|
|
// Initialize of bdf solver. |
129 |
|
✗ |
m_bdfFS.initialize(m_linearProblem[0]->solution()); |
130 |
|
|
// Compute the extrapolate. |
131 |
|
✗ |
m_bdfFS.extrapolate(m_solExtrapolate); |
132 |
|
|
// Initialize it. |
133 |
|
✗ |
m_linearProblem[0]->initExtrapol(m_solExtrapolate); |
134 |
|
|
} else { |
135 |
|
|
// Update bdf. |
136 |
|
✗ |
m_bdfFS.update(m_linearProblem[0]->solution()); |
137 |
|
|
// Compute the extrapolate. |
138 |
|
✗ |
m_bdfFS.extrapolate(m_solExtrapolate); |
139 |
|
|
// Initialize it. |
140 |
|
✗ |
m_linearProblem[0]->updateExtrapol(m_solExtrapolate); |
141 |
|
|
} |
142 |
|
|
} |
143 |
|
|
|
144 |
|
|
// initialize BDF for diffusion |
145 |
|
✗ |
if (iLinPb==1 ) { |
146 |
|
✗ |
if ( m_fstransient->iteration == 1 ) { |
147 |
|
|
// Initialize of bdf solver. |
148 |
|
✗ |
m_bdfFS.initialize(m_linearProblem[1]->solution()); |
149 |
|
|
// Compute the extrapolate. |
150 |
|
✗ |
m_bdfFS.extrapolate(m_solExtrapolate); |
151 |
|
|
// Initialize it. |
152 |
|
✗ |
m_linearProblem[1]->initExtrapol(m_solExtrapolate); |
153 |
|
|
} else { |
154 |
|
|
// Update bdf. |
155 |
|
✗ |
m_bdfFS.update(m_linearProblem[1]->solution()); |
156 |
|
|
// Compute the extrapolate. |
157 |
|
✗ |
m_bdfFS.extrapolate(m_solExtrapolate); |
158 |
|
|
// Initialize it. |
159 |
|
✗ |
m_linearProblem[1]->updateExtrapol(m_solExtrapolate); |
160 |
|
|
} |
161 |
|
|
} |
162 |
|
|
|
163 |
|
|
// last linear problem (projection) |
164 |
|
|
// lumpedModelBC bc computation |
165 |
|
✗ |
if (m_lumpedModelBC && iLinPb==2) { |
166 |
|
✗ |
if ( m_fstransient->iteration == 1 ) { |
167 |
|
|
// m_labelFluxListLumpedModelBC computation with LumpedModelBC labels |
168 |
|
✗ |
m_labelFluxListLumpedModelBC.clear(); |
169 |
|
✗ |
for(std::size_t i=0; i<m_lumpedModelBC->size(); i++) { |
170 |
|
✗ |
m_labelFluxListLumpedModelBC.push_back(FelisceParam::instance().lumpedModelBCLabel[i]); |
171 |
|
|
} |
172 |
|
|
} |
173 |
|
|
// update the fluxes |
174 |
|
|
// justine's warning (10-09-12) : with a velocity with divergence not necesary free... ? |
175 |
|
✗ |
m_linearProblem[0]->computeMeanQuantity(velocity,m_labelFluxListLumpedModelBC,m_fluxLumpedModelBC); |
176 |
|
✗ |
for(std::size_t i=0; i<m_lumpedModelBC->size(); i++) { |
177 |
|
✗ |
m_lumpedModelBC->Q(i) = m_fluxLumpedModelBC[i] ; |
178 |
|
|
} |
179 |
|
|
// iterate the lumpedModelBC parameters |
180 |
|
✗ |
m_lumpedModelBC->iterate(); |
181 |
|
|
} |
182 |
|
|
// Compute the RHS linked to the time scheme (bdf). |
183 |
|
✗ |
std::cout << " ... computeRHSTime - problem " << iLinPb << std::endl; |
184 |
|
✗ |
m_bdfFS.computeRHSTime(m_fstransient->timeStep); |
185 |
|
|
// Fill m_seqBdfRHS. |
186 |
|
✗ |
std::cout << " ... gatherVectorBeforeAssembleMatrixRHS - problem " << iLinPb << std::endl; |
187 |
|
✗ |
m_linearProblem[iLinPb]->gatherVectorBeforeAssembleMatrixRHS(); |
188 |
|
|
//Assembly loop on elements. |
189 |
|
✗ |
std::cout << " ... assembleMatrixRHS() - problem " << iLinPb << std::endl; |
190 |
|
✗ |
m_linearProblem[iLinPb]->assembleMatrixRHS(MpiInfo::rankProc()); |
191 |
|
✗ |
std::cout << " end assembling " << std::endl; |
192 |
|
|
|
193 |
|
|
//Specific operations before solve the system. |
194 |
|
✗ |
postAssemblingMatrixRHS(iLinPb); |
195 |
|
|
//Apply essential transient boundary conditions. |
196 |
|
✗ |
m_linearProblem[iLinPb]->finalizeEssBCTransient(); |
197 |
|
|
// Apply all the boundary conditions. |
198 |
|
✗ |
m_linearProblem[iLinPb]->applyBC(FelisceParam::instance().essentialBoundaryConditionsMethod, MpiInfo::rankProc()); |
199 |
|
✗ |
if (iLinPb==0) { |
200 |
|
✗ |
if (m_fstransient->iteration==2) { |
201 |
|
✗ |
std::cout << " ipb = " << iLinPb << std::endl; |
202 |
|
✗ |
m_linearProblem[iLinPb]->writeMatrixForMatlab(21); |
203 |
|
✗ |
m_linearProblem[iLinPb]->writeRHSForMatlab(21); |
204 |
|
|
//exit(1); |
205 |
|
|
} |
206 |
|
|
} |
207 |
|
|
//Solve linear system. |
208 |
|
✗ |
std::cout << " [*] solve linearProblem " << iLinPb << "[*]" << std::endl; |
209 |
|
✗ |
m_linearProblem[iLinPb]->solve(MpiInfo::rankProc(), MpiInfo::numProc()); |
210 |
|
✗ |
m_linearProblem[iLinPb]->gatherSolution(); |
211 |
|
|
|
212 |
|
|
} |
213 |
|
|
} |
214 |
|
|
|
215 |
|
✗ |
int NSFracStepExplicitModel::getNstate() const { |
216 |
|
✗ |
return m_linearProblem[0]->numDof(); |
217 |
|
|
} |
218 |
|
|
|
219 |
|
✗ |
void NSFracStepExplicitModel::getState(double* & state) { |
220 |
|
✗ |
m_linearProblem[0]->getSolution(state, MpiInfo::numProc(), MpiInfo::rankProc()); |
221 |
|
|
} |
222 |
|
|
|
223 |
|
✗ |
felInt NSFracStepExplicitModel::numDof_swig() { |
224 |
|
✗ |
return m_linearProblem[0]->numDof(); |
225 |
|
|
} |
226 |
|
|
|
227 |
|
✗ |
void NSFracStepExplicitModel::getState_swig(double* data, felInt numDof) { |
228 |
|
|
double * state; |
229 |
|
✗ |
m_linearProblem[0]->getSolution(state, MpiInfo::numProc(), MpiInfo::rankProc()); |
230 |
|
✗ |
for (felInt i=0 ; i<numDof ; i++) { |
231 |
|
✗ |
data[i] = state[i]; |
232 |
|
|
} |
233 |
|
|
} |
234 |
|
|
|
235 |
|
✗ |
void NSFracStepExplicitModel::setState(double* & state) { |
236 |
|
✗ |
m_linearProblem[0]->setSolution(state, MpiInfo::numProc(), MpiInfo::rankProc()); |
237 |
|
|
} |
238 |
|
|
|
239 |
|
|
} |
240 |
|
|
|