Directory: | ./ |
---|---|
File: | PETScInterface/petscVector_inline.hpp |
Date: | 2024-04-14 07:32:34 |
Exec | Total | Coverage | |
---|---|---|---|
Lines: | 332 | 368 | 90.2% |
Branches: | 128 | 406 | 31.5% |
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 | ||
13 | // Dominique Chapelle | ||
14 | // Miguel Fernandez | ||
15 | // Vicente Mataix Ferrandiz | ||
16 | // | ||
17 | |||
18 | #ifndef PETSCVECTOR_HPP | ||
19 | #error Include petscVector.hpp instead of petscVector_inline.hpp | ||
20 | #endif | ||
21 | |||
22 | // System includes | ||
23 | #include <cstdlib> | ||
24 | |||
25 | // External includes | ||
26 | |||
27 | // Project includes | ||
28 | #include "Core/felisce.hpp" | ||
29 | #include "Core/felisceParam.hpp" | ||
30 | #include "PETScInterface/petscMatrix.hpp" | ||
31 | |||
32 | namespace felisce | ||
33 | { | ||
34 | |||
35 | //======================= | ||
36 | // CONSTRUCTOR/DESTRUCTOR | ||
37 | //======================= | ||
38 | 14379 | inline PetscVector::PetscVector() | |
39 | { | ||
40 | 14379 | m_self = nullptr; | |
41 | 14379 | } | |
42 | |||
43 | /***********************************************************************************/ | ||
44 | /***********************************************************************************/ | ||
45 | |||
46 | 120 | inline PetscVector::PetscVector(const PetscVector& b) | |
47 | { | ||
48 |
2/2✓ Branch 1 taken 104 times.
✓ Branch 2 taken 16 times.
|
120 | if ( b.isNotNull() ) { |
49 | 104 | duplicateFrom(b); | |
50 | 104 | copyFrom(b); | |
51 | } | ||
52 | 120 | } | |
53 | |||
54 | /***********************************************************************************/ | ||
55 | /***********************************************************************************/ | ||
56 | |||
57 | 195 | inline PetscVector& PetscVector::operator=(const PetscVector& b) | |
58 | { | ||
59 | 195 | this->destroy(); | |
60 | 195 | PetscObjectReference((PetscObject)b.m_self); | |
61 | 195 | m_self = b.m_self; | |
62 | 195 | return *this; | |
63 | } | ||
64 | |||
65 | /***********************************************************************************/ | ||
66 | /***********************************************************************************/ | ||
67 | |||
68 | 4 | inline PetscVector PetscVector::operator*(const PetscInt scalar) | |
69 | { | ||
70 | 4 | PetscVector a(*this); | |
71 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(PetscScalar(scalar)); |
72 | 4 | return a; | |
73 | } | ||
74 | /***********************************************************************************/ | ||
75 | /***********************************************************************************/ | ||
76 | |||
77 | 4 | inline PetscVector PetscVector::operator*(const PetscInt scalar) const | |
78 | { | ||
79 | 4 | PetscVector a(*this); | |
80 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(PetscScalar(scalar)); |
81 | 4 | return a; | |
82 | } | ||
83 | |||
84 | /***********************************************************************************/ | ||
85 | /***********************************************************************************/ | ||
86 | |||
87 | 4 | inline PetscVector PetscVector::operator*(const PetscScalar scalar) | |
88 | { | ||
89 | 4 | PetscVector a(*this); | |
90 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(scalar); |
91 | 4 | return a; | |
92 | } | ||
93 | /***********************************************************************************/ | ||
94 | /***********************************************************************************/ | ||
95 | |||
96 | 4 | inline PetscVector PetscVector::operator*(const PetscScalar scalar) const | |
97 | { | ||
98 | 4 | PetscVector a(*this); | |
99 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(scalar); |
100 | 4 | return a; | |
101 | } | ||
102 | |||
103 | /***********************************************************************************/ | ||
104 | /***********************************************************************************/ | ||
105 | |||
106 | 8 | inline PetscScalar PetscVector::operator*(const PetscVector& b) | |
107 | { | ||
108 | PetscScalar a; | ||
109 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | dot(*this, b, &a); |
110 | 8 | return a; | |
111 | } | ||
112 | |||
113 | /***********************************************************************************/ | ||
114 | /***********************************************************************************/ | ||
115 | |||
116 | 8 | inline PetscScalar PetscVector::operator*(const PetscVector& b) const | |
117 | { | ||
118 | PetscScalar a; | ||
119 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | dot(*this, b, &a); |
120 | 8 | return a; | |
121 | } | ||
122 | |||
123 | /***********************************************************************************/ | ||
124 | /***********************************************************************************/ | ||
125 | |||
126 | 4 | inline PetscVector& PetscVector::operator*=(const PetscInt scalar) | |
127 | { | ||
128 | 4 | scale(PetscScalar(scalar)); | |
129 | 4 | return *this; | |
130 | } | ||
131 | |||
132 | /***********************************************************************************/ | ||
133 | /***********************************************************************************/ | ||
134 | |||
135 | 4 | inline PetscVector& PetscVector::operator*=(const PetscScalar scalar) | |
136 | { | ||
137 | 4 | scale(scalar); | |
138 | 4 | return *this; | |
139 | } | ||
140 | |||
141 | /***********************************************************************************/ | ||
142 | /***********************************************************************************/ | ||
143 | |||
144 | 4 | inline PetscVector PetscVector::operator/(const PetscInt scalar) | |
145 | { | ||
146 | 4 | PetscVector a(*this); | |
147 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(1.0/PetscScalar(scalar)); |
148 | 4 | return a; | |
149 | } | ||
150 | |||
151 | /***********************************************************************************/ | ||
152 | /***********************************************************************************/ | ||
153 | |||
154 | 4 | inline PetscVector PetscVector::operator/(const PetscInt scalar) const | |
155 | { | ||
156 | 4 | PetscVector a(*this); | |
157 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(1.0/PetscScalar(scalar)); |
158 | 4 | return a; | |
159 | } | ||
160 | |||
161 | /***********************************************************************************/ | ||
162 | /***********************************************************************************/ | ||
163 | |||
164 | 4 | inline PetscVector PetscVector::operator/(const PetscScalar scalar) | |
165 | { | ||
166 | 4 | PetscVector a(*this); | |
167 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(1.0/scalar); |
168 | 4 | return a; | |
169 | } | ||
170 | |||
171 | /***********************************************************************************/ | ||
172 | /***********************************************************************************/ | ||
173 | |||
174 | 4 | inline PetscVector PetscVector::operator/(const PetscScalar scalar) const | |
175 | { | ||
176 | 4 | PetscVector a(*this); | |
177 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.scale(1.0/scalar); |
178 | 4 | return a; | |
179 | } | ||
180 | |||
181 | /***********************************************************************************/ | ||
182 | /***********************************************************************************/ | ||
183 | |||
184 | 4 | inline PetscVector& PetscVector::operator/=(const PetscInt scalar) | |
185 | { | ||
186 | 4 | scale(1.0/PetscScalar(scalar)); | |
187 | 4 | return *this; | |
188 | } | ||
189 | |||
190 | /***********************************************************************************/ | ||
191 | /***********************************************************************************/ | ||
192 | |||
193 | 4 | inline PetscVector& PetscVector::operator/=(const PetscScalar scalar) | |
194 | { | ||
195 | 4 | scale(1.0/scalar); | |
196 | 4 | return *this; | |
197 | } | ||
198 | |||
199 | /***********************************************************************************/ | ||
200 | /***********************************************************************************/ | ||
201 | |||
202 | 8 | inline PetscVector PetscVector::operator+(const PetscVector& b) | |
203 | { | ||
204 | 8 | PetscVector a(*this); | |
205 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | a.axpy(b); |
206 | 8 | return a; | |
207 | } | ||
208 | |||
209 | /***********************************************************************************/ | ||
210 | /***********************************************************************************/ | ||
211 | |||
212 | 4 | inline PetscVector PetscVector::operator+(const PetscVector& b) const | |
213 | { | ||
214 | 4 | PetscVector a(*this); | |
215 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.axpy(b); |
216 | 4 | return a; | |
217 | } | ||
218 | |||
219 | /***********************************************************************************/ | ||
220 | /***********************************************************************************/ | ||
221 | |||
222 | 4 | inline PetscVector& PetscVector::operator+=(const PetscVector& b) | |
223 | { | ||
224 | 4 | axpy(b); | |
225 | 4 | return *this; | |
226 | } | ||
227 | |||
228 | /***********************************************************************************/ | ||
229 | /***********************************************************************************/ | ||
230 | |||
231 | 4 | inline PetscVector PetscVector::operator-(const PetscVector& b) | |
232 | { | ||
233 | 4 | PetscVector a(*this); | |
234 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.axpy(b, -1.0); |
235 | 4 | return a; | |
236 | } | ||
237 | |||
238 | /***********************************************************************************/ | ||
239 | /***********************************************************************************/ | ||
240 | |||
241 | 4 | inline PetscVector PetscVector::operator-(const PetscVector& b) const | |
242 | { | ||
243 | 4 | PetscVector a(*this); | |
244 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | a.axpby(b, -1.0); |
245 | 4 | return a; | |
246 | } | ||
247 | |||
248 | /***********************************************************************************/ | ||
249 | /***********************************************************************************/ | ||
250 | |||
251 | 4 | inline PetscVector& PetscVector::operator-=(const PetscVector& b) | |
252 | { | ||
253 | 4 | axpy(b, -1.0); | |
254 | 4 | return *this; | |
255 | } | ||
256 | |||
257 | /***********************************************************************************/ | ||
258 | /***********************************************************************************/ | ||
259 | |||
260 | 202 | inline PetscScalar PetscVector::operator[] (const PetscInt Index) | |
261 | { | ||
262 | double value; | ||
263 |
1/2✓ Branch 1 taken 202 times.
✗ Branch 2 not taken.
|
202 | PetscErrorCode code = getValue(Index, value); |
264 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 202 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
202 | CHKERRQ(code); |
265 | 202 | return value; | |
266 | } | ||
267 | |||
268 | /***********************************************************************************/ | ||
269 | /***********************************************************************************/ | ||
270 | |||
271 | 40 | inline PetscScalar PetscVector::operator[] (const PetscInt Index) const | |
272 | { | ||
273 | double value; | ||
274 |
1/2✓ Branch 1 taken 40 times.
✗ Branch 2 not taken.
|
40 | PetscErrorCode code = getValue(Index, value); |
275 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 40 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
40 | CHKERRQ(code); |
276 | 40 | return value; | |
277 | } | ||
278 | |||
279 | /***********************************************************************************/ | ||
280 | /***********************************************************************************/ | ||
281 | |||
282 | 14499 | inline PetscVector::~PetscVector() | |
283 | { | ||
284 | // In PETSc, Vec are handled like shared pointers. | ||
285 | 14499 | this->destroy(); | |
286 | 14499 | } | |
287 | |||
288 | //================= | ||
289 | // STATIC FUNCTIONS | ||
290 | //================= | ||
291 | |||
292 | 138 | inline PetscVector PetscVector::null() | |
293 | { | ||
294 | 138 | PetscVector nullVector; | |
295 | 138 | return nullVector; | |
296 | } | ||
297 | |||
298 | //=================== | ||
299 | // CREATION FUNCTIONS | ||
300 | //=================== | ||
301 | |||
302 | 187 | inline PetscErrorCode PetscVector::create(MPI_Comm comm) | |
303 | { | ||
304 | |||
305 | 187 | PetscErrorCode code = VecCreate(comm,&m_self); | |
306 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 187 times.
|
187 | CHKERRQ(code); |
307 | 187 | return code; | |
308 | } | ||
309 | |||
310 | /***********************************************************************************/ | ||
311 | /***********************************************************************************/ | ||
312 | |||
313 | 857 | inline PetscErrorCode PetscVector::createSeq(MPI_Comm comm,PetscInt n) | |
314 | { | ||
315 | 857 | PetscErrorCode code = VecCreateSeq(comm,n,&m_self); | |
316 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 857 times.
|
857 | CHKERRQ(code); |
317 | 857 | return code; | |
318 | } | ||
319 | |||
320 | /***********************************************************************************/ | ||
321 | /***********************************************************************************/ | ||
322 | |||
323 | 1974 | inline PetscErrorCode PetscVector::createMPI(MPI_Comm comm,PetscInt n,PetscInt N) | |
324 | { | ||
325 | 1974 | PetscErrorCode code = VecCreateMPI(comm,n,N,&m_self); | |
326 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1974 times.
|
1974 | CHKERRQ(code); |
327 | 1974 | return code; | |
328 | } | ||
329 | |||
330 | //====================== | ||
331 | // DESTRUCTION FUNCTIONS | ||
332 | //====================== | ||
333 | |||
334 | 17204 | inline PetscErrorCode PetscVector::destroy() | |
335 | { | ||
336 | 17204 | PetscErrorCode code = VecDestroy(&m_self); | |
337 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 17204 times.
|
17204 | CHKERRQ(code); |
338 | 17204 | return code; | |
339 | } | ||
340 | |||
341 | //=============== | ||
342 | // COPY FUNCTIONS | ||
343 | //=============== | ||
344 | |||
345 | 38164 | inline PetscErrorCode PetscVector::copyFrom(const PetscVector& vector) const | |
346 | { | ||
347 | 38164 | const PetscErrorCode code = VecCopy(vector.m_self,m_self); | |
348 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 38164 times.
|
38164 | CHKERRQ(code); |
349 | 38164 | return code; | |
350 | } | ||
351 | |||
352 | /***********************************************************************************/ | ||
353 | /***********************************************************************************/ | ||
354 | |||
355 | 6073 | inline PetscErrorCode PetscVector::duplicateFrom(const PetscVector& vector) | |
356 | { | ||
357 |
2/2✓ Branch 1 taken 150 times.
✓ Branch 2 taken 5923 times.
|
6073 | if ( this->isNotNull() ) |
358 | 150 | this->destroy(); | |
359 | 6073 | const PetscErrorCode code = VecDuplicate(vector.m_self,&m_self); | |
360 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6073 times.
|
6073 | CHKERRQ(code); |
361 | 6073 | return code; | |
362 | } | ||
363 | |||
364 | //========================= | ||
365 | // MODIFICATION FUNCTIONS | ||
366 | //========================= | ||
367 | |||
368 | 92508 | inline PetscErrorCode PetscVector::assembly() const | |
369 | { | ||
370 | 92508 | PetscErrorCode code = VecAssemblyBegin(m_self); | |
371 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 92508 times.
|
92508 | CHKERRQ(code); |
372 | 92508 | code = VecAssemblyEnd(m_self); | |
373 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 92508 times.
|
92508 | CHKERRQ(code); |
374 | 92508 | return code; | |
375 | } | ||
376 | |||
377 | /***********************************************************************************/ | ||
378 | /***********************************************************************************/ | ||
379 | |||
380 | 6312700 | inline PetscErrorCode PetscVector::restoreArray(PetscScalar **a) const | |
381 | { | ||
382 | 6312700 | const PetscErrorCode code = VecRestoreArray(m_self,a); | |
383 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6312700 times.
|
6312700 | CHKERRQ(code); |
384 | 6312700 | return code; | |
385 | } | ||
386 | |||
387 | /***********************************************************************************/ | ||
388 | /***********************************************************************************/ | ||
389 | |||
390 | 3261 | inline PetscErrorCode PetscVector::set(PetscScalar const & value) const | |
391 | { | ||
392 | 3261 | const PetscErrorCode code = VecSet(m_self,value); | |
393 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 3261 times.
|
3261 | CHKERRQ(code); |
394 | 3261 | return code; | |
395 | } | ||
396 | |||
397 | /***********************************************************************************/ | ||
398 | /***********************************************************************************/ | ||
399 | |||
400 | 349 | inline PetscErrorCode PetscVector::setFromOptions() const | |
401 | { | ||
402 | 349 | const PetscErrorCode code = VecSetFromOptions(m_self); | |
403 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 349 times.
|
349 | CHKERRQ(code); |
404 | 349 | return code; | |
405 | } | ||
406 | |||
407 | /***********************************************************************************/ | ||
408 | /***********************************************************************************/ | ||
409 | |||
410 | inline PetscErrorCode PetscVector::setOption(VecOption op, PetscBool flg) const | ||
411 | { | ||
412 | const PetscErrorCode code = VecSetOption(m_self,op,flg); | ||
413 | CHKERRQ(code); | ||
414 | return code; | ||
415 | } | ||
416 | |||
417 | /***********************************************************************************/ | ||
418 | /***********************************************************************************/ | ||
419 | |||
420 | 187 | inline PetscErrorCode PetscVector::setSizes(PetscInt n,PetscInt N) const | |
421 | { | ||
422 | 187 | const PetscErrorCode code = VecSetSizes(m_self,n,N); | |
423 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 187 times.
|
187 | CHKERRQ(code); |
424 | 187 | return code; | |
425 | } | ||
426 | |||
427 | /***********************************************************************************/ | ||
428 | /***********************************************************************************/ | ||
429 | |||
430 | 150 | inline PetscErrorCode PetscVector::setType(VecType method) const | |
431 | { | ||
432 | 150 | const PetscErrorCode code = VecSetType(m_self,method); | |
433 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 150 times.
|
150 | CHKERRQ(code); |
434 | 150 | return code; | |
435 | } | ||
436 | |||
437 | /***********************************************************************************/ | ||
438 | /***********************************************************************************/ | ||
439 | |||
440 | 6059113 | inline PetscErrorCode PetscVector::setValue(PetscInt row,PetscScalar value, InsertMode mode) const | |
441 | { | ||
442 | 6059113 | const PetscErrorCode code = VecSetValue(m_self,row,value,mode); | |
443 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6059113 times.
|
6059113 | CHKERRQ(code); |
444 | 6059113 | return code; | |
445 | } | ||
446 | |||
447 | /***********************************************************************************/ | ||
448 | /***********************************************************************************/ | ||
449 | |||
450 | 3064 | inline PetscErrorCode PetscVector::setValueOnce(PetscInt row,PetscScalar value, InsertMode mode) const | |
451 | { | ||
452 | 3064 | PetscErrorCode code = 0; | |
453 | // Get rank | ||
454 | int rank; | ||
455 |
1/2✓ Branch 1 taken 3064 times.
✗ Branch 2 not taken.
|
3064 | MPI_Comm_rank(PETSC_COMM_WORLD, &rank); |
456 | // Set value | ||
457 |
2/2✓ Branch 0 taken 790 times.
✓ Branch 1 taken 2274 times.
|
3064 | if (rank == 0) { |
458 |
1/2✓ Branch 1 taken 790 times.
✗ Branch 2 not taken.
|
790 | code = VecSetValue(m_self,row,value,mode); |
459 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 790 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
790 | CHKERRQ(code); |
460 | } | ||
461 | 3064 | return code; | |
462 | } | ||
463 | |||
464 | /***********************************************************************************/ | ||
465 | /***********************************************************************************/ | ||
466 | |||
467 | 23388674 | inline PetscErrorCode PetscVector::setValues(PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora) const | |
468 | { | ||
469 | 23388674 | const PetscErrorCode code = VecSetValues(m_self,ni,ix,y,iora); | |
470 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 23388674 times.
|
23388674 | CHKERRQ(code); |
471 | 23388674 | return code; | |
472 | } | ||
473 | |||
474 | /***********************************************************************************/ | ||
475 | /***********************************************************************************/ | ||
476 | |||
477 | 4 | inline PetscErrorCode PetscVector::setValuesOnce(PetscInt ni,const PetscInt ix[],const PetscScalar y[],InsertMode iora) const | |
478 | { | ||
479 | 4 | PetscErrorCode code = 0; | |
480 | // Get rank | ||
481 | int rank; | ||
482 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | MPI_Comm_rank(PETSC_COMM_WORLD, &rank); |
483 | // Set value | ||
484 |
2/2✓ Branch 0 taken 1 times.
✓ Branch 1 taken 3 times.
|
4 | if (rank == 0) { |
485 |
1/2✓ Branch 1 taken 1 times.
✗ Branch 2 not taken.
|
1 | code = VecSetValues(m_self,ni,ix,y,iora); |
486 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 1 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
1 | CHKERRQ(code); |
487 | } | ||
488 | 4 | return code; | |
489 | } | ||
490 | |||
491 | /***********************************************************************************/ | ||
492 | /***********************************************************************************/ | ||
493 | |||
494 | 79863 | inline PetscErrorCode PetscVector::zeroEntries() const | |
495 | { | ||
496 | 79863 | const PetscErrorCode code = VecZeroEntries(m_self); | |
497 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 79863 times.
|
79863 | CHKERRQ(code); |
498 | 79863 | return code; | |
499 | } | ||
500 | |||
501 | //================= | ||
502 | // ACCESS FUNCTIONS | ||
503 | //================= | ||
504 | |||
505 | 6312704 | inline PetscErrorCode PetscVector::getArray(PetscScalar** array) const | |
506 | { | ||
507 | 6312704 | const PetscErrorCode code = VecGetArray(m_self,array); | |
508 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6312704 times.
|
6312704 | CHKERRQ(code); |
509 | 6312704 | return code; | |
510 | } | ||
511 | |||
512 | /***********************************************************************************/ | ||
513 | /***********************************************************************************/ | ||
514 | |||
515 | 4 | inline PetscScalar* PetscVector::getArray() const | |
516 | { | ||
517 | PetscScalar* array; | ||
518 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | getArray(&array); |
519 | 4 | return array; | |
520 | } | ||
521 | |||
522 | /***********************************************************************************/ | ||
523 | /***********************************************************************************/ | ||
524 | |||
525 | 4425 | inline PetscErrorCode PetscVector::getAllValuesInAppOrdering(const AO& ao, std::vector<double>& array) const | |
526 | { | ||
527 | PetscErrorCode code; | ||
528 | felInt numDof; | ||
529 | 8850 | code = this->getSize(&numDof); | |
530 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 4425 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
4425 | CHKERRQ(code); |
531 |
1/2✓ Branch 1 taken 4425 times.
✗ Branch 2 not taken.
|
4425 | array.resize(numDof); |
532 | 4425 | felInt identityMap[numDof]; | |
533 |
2/2✓ Branch 0 taken 16273911 times.
✓ Branch 1 taken 4425 times.
|
16278336 | for (felInt i=0; i<numDof; i++) |
534 | 16273911 | identityMap[i] = i; | |
535 |
1/2✓ Branch 1 taken 4425 times.
✗ Branch 2 not taken.
|
4425 | AOApplicationToPetsc(ao,numDof,identityMap); |
536 |
1/2✓ Branch 2 taken 4425 times.
✗ Branch 3 not taken.
|
4425 | code = this->getValues(numDof,identityMap,array.data()); |
537 | 4425 | return code; | |
538 |
1/2✓ Branch 1 taken 4425 times.
✗ Branch 2 not taken.
|
4425 | } |
539 | |||
540 | /***********************************************************************************/ | ||
541 | /***********************************************************************************/ | ||
542 | |||
543 | 6312653 | inline PetscErrorCode PetscVector::getLocalSize(PetscInt* localSize) const | |
544 | { | ||
545 | 6312653 | const PetscErrorCode code = VecGetLocalSize(m_self,localSize); | |
546 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 6312653 times.
|
6312653 | CHKERRQ(code); |
547 | 6312653 | return code; | |
548 | } | ||
549 | |||
550 | /***********************************************************************************/ | ||
551 | /***********************************************************************************/ | ||
552 | |||
553 | 8 | inline PetscInt PetscVector::getLocalSize() const | |
554 | { | ||
555 | PetscInt local_size; | ||
556 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | getLocalSize(&local_size); |
557 | 8 | return local_size; | |
558 | } | ||
559 | |||
560 | /***********************************************************************************/ | ||
561 | /***********************************************************************************/ | ||
562 | |||
563 | 4966 | inline PetscErrorCode PetscVector::getSize(PetscInt* size) const | |
564 | { | ||
565 | 4966 | const PetscErrorCode code = VecGetSize(m_self,size); | |
566 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 4966 times.
|
4966 | CHKERRQ(code); |
567 | 4966 | return code; | |
568 | } | ||
569 | |||
570 | /***********************************************************************************/ | ||
571 | /***********************************************************************************/ | ||
572 | |||
573 | 32 | inline PetscInt PetscVector::getSize() const | |
574 | { | ||
575 | PetscInt size; | ||
576 |
1/2✓ Branch 1 taken 32 times.
✗ Branch 2 not taken.
|
32 | getSize(&size); |
577 | 32 | return size; | |
578 | } | ||
579 | |||
580 | /***********************************************************************************/ | ||
581 | /***********************************************************************************/ | ||
582 | |||
583 | 69055536 | inline PetscErrorCode PetscVector::getValues(PetscInt ni,const PetscInt ix[],PetscScalar y[]) const | |
584 | { | ||
585 | 69055536 | const PetscErrorCode code = VecGetValues(m_self,ni,ix,y); | |
586 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 69055536 times.
|
69055536 | CHKERRQ(code); |
587 | 69055536 | return code; | |
588 | } | ||
589 | |||
590 | /***********************************************************************************/ | ||
591 | /***********************************************************************************/ | ||
592 | |||
593 | 262 | inline PetscErrorCode PetscVector::getValue(const PetscInt ix,PetscScalar& y) const | |
594 | { | ||
595 | 262 | return getValues(1, &ix, &y); | |
596 | } | ||
597 | |||
598 | //================== | ||
599 | // BOOLEAN FUNCTIONS | ||
600 | //================== | ||
601 | |||
602 | 513 | inline bool PetscVector::isNull() | |
603 | { | ||
604 | 513 | return m_self == nullptr; | |
605 | } | ||
606 | |||
607 | /***********************************************************************************/ | ||
608 | /***********************************************************************************/ | ||
609 | |||
610 | inline bool PetscVector::isNull() const | ||
611 | { | ||
612 | return m_self == nullptr; | ||
613 | } | ||
614 | |||
615 | /***********************************************************************************/ | ||
616 | /***********************************************************************************/ | ||
617 | |||
618 | 103207 | inline bool PetscVector::isNotNull() | |
619 | { | ||
620 | 103207 | return m_self != nullptr; | |
621 | } | ||
622 | |||
623 | /***********************************************************************************/ | ||
624 | /***********************************************************************************/ | ||
625 | |||
626 | 156 | inline bool PetscVector::isNotNull() const | |
627 | { | ||
628 | 156 | return m_self != nullptr; | |
629 | } | ||
630 | |||
631 | //======================= | ||
632 | // MATHEMATICAL FUNCTIONS | ||
633 | //======================= | ||
634 | ✗ | inline PetscErrorCode PetscVector::axpbypcz(PetscScalar alpha, PetscScalar beta, PetscScalar gamma, const PetscVector& x, const PetscVector& y) const | |
635 | { | ||
636 | ✗ | const PetscErrorCode code = VecAXPBYPCZ(m_self,alpha,beta,gamma,x.m_self, y.m_self); | |
637 | ✗ | CHKERRQ(code); | |
638 | ✗ | return code; | |
639 | } | ||
640 | |||
641 | /***********************************************************************************/ | ||
642 | /***********************************************************************************/ | ||
643 | |||
644 | 15190 | inline PetscErrorCode PetscVector::axpby(PetscScalar alpha, PetscScalar beta, const PetscVector& x) const | |
645 | { | ||
646 | 15190 | const PetscErrorCode code = VecAXPBY(m_self,alpha,beta,x.m_self); | |
647 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 15190 times.
|
15190 | CHKERRQ(code); |
648 | 15190 | return code; | |
649 | } | ||
650 | |||
651 | /***********************************************************************************/ | ||
652 | /***********************************************************************************/ | ||
653 | |||
654 | 41894 | inline PetscErrorCode PetscVector::axpy(PetscScalar alpha, const PetscVector& x) const | |
655 | { | ||
656 | // m_self = m_self + alpha * x.m_self | ||
657 | 41894 | const PetscErrorCode code = VecAXPY(m_self,alpha,x.m_self); | |
658 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 41894 times.
|
41894 | CHKERRQ(code); |
659 | 41894 | return code; | |
660 | } | ||
661 | |||
662 | /***********************************************************************************/ | ||
663 | /***********************************************************************************/ | ||
664 | |||
665 | 41901 | inline PetscErrorCode PetscVector::max(PetscReal * val, PetscInt* p) const | |
666 | { | ||
667 | 41901 | const PetscErrorCode code = VecMax(m_self,p,val); | |
668 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 41901 times.
|
41901 | CHKERRQ(code); |
669 | 41901 | return code; | |
670 | } | ||
671 | |||
672 | /***********************************************************************************/ | ||
673 | /***********************************************************************************/ | ||
674 | |||
675 | inline PetscErrorCode PetscVector::max(PetscInt * p,PetscReal * val) const | ||
676 | { | ||
677 | const PetscErrorCode code = VecMax(m_self,p,val); | ||
678 | CHKERRQ(code); | ||
679 | return code; | ||
680 | } | ||
681 | /***********************************************************************************/ | ||
682 | /***********************************************************************************/ | ||
683 | |||
684 | 4 | inline PetscReal PetscVector::max() const | |
685 | { | ||
686 | PetscReal max_value; | ||
687 |
1/2✓ Branch 1 taken 4 times.
✗ Branch 2 not taken.
|
4 | max(&max_value); |
688 | 4 | return max_value; | |
689 | } | ||
690 | |||
691 | /***********************************************************************************/ | ||
692 | /***********************************************************************************/ | ||
693 | |||
694 | 41897 | inline PetscErrorCode PetscVector::min(PetscReal * val, PetscInt* p) const | |
695 | { | ||
696 | 41897 | const PetscErrorCode code = VecMin(m_self,p,val); | |
697 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 41897 times.
|
41897 | CHKERRQ(code); |
698 | 41897 | return code; | |
699 | } | ||
700 | |||
701 | /***********************************************************************************/ | ||
702 | /***********************************************************************************/ | ||
703 | |||
704 | inline PetscErrorCode PetscVector::min(PetscInt * p,PetscReal * val) const | ||
705 | { | ||
706 | const PetscErrorCode code = VecMin(m_self,p,val); | ||
707 | CHKERRQ(code); | ||
708 | return code; | ||
709 | } | ||
710 | |||
711 | /***********************************************************************************/ | ||
712 | /***********************************************************************************/ | ||
713 | |||
714 | inline PetscReal PetscVector::min() const | ||
715 | { | ||
716 | PetscReal min_value; | ||
717 | min(&min_value); | ||
718 | return min_value; | ||
719 | } | ||
720 | |||
721 | /***********************************************************************************/ | ||
722 | /***********************************************************************************/ | ||
723 | |||
724 | 33 | inline PetscErrorCode PetscVector::norm(NormType type, PetscReal* val) const | |
725 | { | ||
726 | 33 | const PetscErrorCode code = VecNorm(m_self,type,val); | |
727 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 33 times.
|
33 | CHKERRQ(code); |
728 | 33 | return code; | |
729 | } | ||
730 | |||
731 | /***********************************************************************************/ | ||
732 | /***********************************************************************************/ | ||
733 | |||
734 | 8 | inline PetscReal PetscVector::norm(const NormType type) const | |
735 | { | ||
736 | 8 | PetscReal norm_value = 0.0; | |
737 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | norm(type, &norm_value); |
738 | 8 | return norm_value; | |
739 | } | ||
740 | |||
741 | /***********************************************************************************/ | ||
742 | /***********************************************************************************/ | ||
743 | |||
744 | 1453 | inline PetscErrorCode PetscVector::scale(PetscScalar scalar) const | |
745 | { | ||
746 | 1453 | const PetscErrorCode code = VecScale(m_self,scalar); | |
747 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1453 times.
|
1453 | CHKERRQ(code); |
748 | 1453 | return code; | |
749 | } | ||
750 | |||
751 | /***********************************************************************************/ | ||
752 | /***********************************************************************************/ | ||
753 | |||
754 | 100 | inline PetscErrorCode PetscVector::shift(PetscScalar scalar) const | |
755 | { | ||
756 | 100 | const PetscErrorCode code = VecShift(m_self,scalar); | |
757 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 100 times.
|
100 | CHKERRQ(code); |
758 | 100 | return code; | |
759 | } | ||
760 | |||
761 | /***********************************************************************************/ | ||
762 | /***********************************************************************************/ | ||
763 | |||
764 | 8 | inline PetscErrorCode PetscVector::sum(PetscScalar *sum) const | |
765 | { | ||
766 | 8 | const PetscErrorCode code = VecSum(m_self,sum); | |
767 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
|
8 | CHKERRQ(code); |
768 | 8 | return code; | |
769 | } | ||
770 | |||
771 | /***********************************************************************************/ | ||
772 | /***********************************************************************************/ | ||
773 | |||
774 | 8 | inline PetscScalar PetscVector::sum() const | |
775 | { | ||
776 | 8 | PetscReal sum_value = 0.0; | |
777 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | sum(&sum_value); |
778 | 8 | return sum_value; | |
779 | } | ||
780 | |||
781 | /***********************************************************************************/ | ||
782 | /***********************************************************************************/ | ||
783 | |||
784 | 79 | inline PetscErrorCode PetscVector::abs() const | |
785 | { | ||
786 | 79 | const PetscErrorCode code = VecAbs(m_self); | |
787 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 79 times.
|
79 | CHKERRQ(code); |
788 | 79 | return code; | |
789 | } | ||
790 | |||
791 | //======================== | ||
792 | // COMMUNICATION FUNCTIONS | ||
793 | //======================== | ||
794 | |||
795 | 126 | inline PetscErrorCode PetscVector::scatterToAll(PetscVector& y, InsertMode addv, ScatterMode mode) const | |
796 | { | ||
797 | VecScatter inctx; | ||
798 |
1/2✓ Branch 1 taken 126 times.
✗ Branch 2 not taken.
|
126 | PetscErrorCode code = VecScatterCreateToAll(m_self,&inctx,&y.m_self); |
799 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 126 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
126 | CHKERRQ(code); |
800 |
1/2✓ Branch 1 taken 126 times.
✗ Branch 2 not taken.
|
126 | code = VecScatterBegin(inctx,m_self,y.m_self,addv,mode); |
801 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 126 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
126 | CHKERRQ(code); |
802 |
1/2✓ Branch 1 taken 126 times.
✗ Branch 2 not taken.
|
126 | code = VecScatterEnd(inctx,m_self,y.m_self,addv,mode); |
803 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 126 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
126 | CHKERRQ(code); |
804 |
1/2✓ Branch 1 taken 126 times.
✗ Branch 2 not taken.
|
126 | code = VecScatterDestroy(&inctx); |
805 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 126 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
126 | CHKERRQ(code); |
806 | 126 | return code; | |
807 | } | ||
808 | |||
809 | /***********************************************************************************/ | ||
810 | /***********************************************************************************/ | ||
811 | |||
812 | 200 | inline PetscErrorCode PetscVector::broadCastSequentialVector(MPI_Comm comm, int master ) const | |
813 | { | ||
814 | double * array; | ||
815 | int size; | ||
816 |
1/2✓ Branch 1 taken 200 times.
✗ Branch 2 not taken.
|
200 | PetscErrorCode code = this->getSize(&size); |
817 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 200 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
200 | CHKERRQ(code); |
818 |
1/2✓ Branch 1 taken 200 times.
✗ Branch 2 not taken.
|
200 | code = this->getArray(&array); |
819 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 200 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
200 | CHKERRQ(code); |
820 |
1/2✓ Branch 1 taken 200 times.
✗ Branch 2 not taken.
|
200 | MPI_Bcast( array, size, MPI_DOUBLE, master, comm); |
821 |
1/2✓ Branch 1 taken 200 times.
✗ Branch 2 not taken.
|
200 | code = this->restoreArray(&array); |
822 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 200 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
200 | CHKERRQ(code); |
823 | 200 | return code; | |
824 | } | ||
825 | |||
826 | /***********************************************************************************/ | ||
827 | /***********************************************************************************/ | ||
828 | |||
829 | 96972 | inline PetscErrorCode PetscVector::scatterToAllNotCreatingVector(PetscVector& y,InsertMode addv,ScatterMode mode) const | |
830 | { | ||
831 | VecScatter inctx; | ||
832 | PetscErrorCode code; | ||
833 |
2/6✓ Branch 1 taken 96972 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 96972 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
96972 | code = VecScatterCreateToAll(m_self,&inctx,nullptr); CHKERRQ(code); |
834 |
2/6✓ Branch 1 taken 96972 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 96972 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
96972 | code = VecScatterBegin(inctx,m_self,y.m_self,addv,mode); CHKERRQ(code); |
835 |
2/6✓ Branch 1 taken 96972 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 96972 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
96972 | code = VecScatterEnd(inctx,m_self,y.m_self,addv,mode); CHKERRQ(code); |
836 |
2/6✓ Branch 1 taken 96972 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 96972 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
96972 | code = VecScatterDestroy(&inctx); CHKERRQ(code); |
837 | 96972 | return code; | |
838 | } | ||
839 | |||
840 | /***********************************************************************************/ | ||
841 | /***********************************************************************************/ | ||
842 | |||
843 | 631 | inline PetscErrorCode PetscVector::scatterToZero(PetscVector& y,InsertMode addv,ScatterMode mode) const | |
844 | { | ||
845 | VecScatter inctx; | ||
846 |
1/2✓ Branch 1 taken 631 times.
✗ Branch 2 not taken.
|
631 | PetscErrorCode code = VecScatterCreateToZero(m_self,&inctx,&y.m_self); |
847 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 631 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
631 | CHKERRQ(code); |
848 |
1/2✓ Branch 1 taken 631 times.
✗ Branch 2 not taken.
|
631 | code = VecScatterBegin(inctx,m_self,y.m_self,addv,mode); |
849 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 631 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
631 | CHKERRQ(code); |
850 |
1/2✓ Branch 1 taken 631 times.
✗ Branch 2 not taken.
|
631 | code = VecScatterEnd(inctx,m_self,y.m_self,addv,mode); |
851 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 631 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
631 | CHKERRQ(code); |
852 |
1/2✓ Branch 1 taken 631 times.
✗ Branch 2 not taken.
|
631 | code = VecScatterDestroy(&inctx); |
853 |
1/4✗ Branch 0 not taken.
✓ Branch 1 taken 631 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
|
631 | CHKERRQ(code); |
854 | 631 | return code; | |
855 | } | ||
856 | |||
857 | /***********************************************************************************/ | ||
858 | /***********************************************************************************/ | ||
859 | |||
860 | 11 | inline PetscErrorCode PetscVector::scatterToZeroNotCreatingVector(PetscVector& y,InsertMode addv,ScatterMode mode) const | |
861 | { | ||
862 | VecScatter inctx; | ||
863 | PetscErrorCode code; | ||
864 |
2/6✓ Branch 1 taken 11 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 11 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
11 | code = VecScatterCreateToZero(m_self,&inctx,nullptr); CHKERRQ(code); |
865 |
2/6✓ Branch 1 taken 11 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 11 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
11 | code = VecScatterBegin(inctx,m_self,y.m_self,addv,mode); CHKERRQ(code); |
866 |
2/6✓ Branch 1 taken 11 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 11 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
11 | code = VecScatterEnd(inctx,m_self,y.m_self,addv,mode); CHKERRQ(code); |
867 |
2/6✓ Branch 1 taken 11 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 11 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
11 | code = VecScatterDestroy(&inctx); CHKERRQ(code); |
868 | 11 | return code; | |
869 | } | ||
870 | |||
871 | //================= | ||
872 | // OUTPUT FUNCTIONS | ||
873 | //================= | ||
874 | |||
875 | 8 | inline PetscErrorCode PetscVector::view(const PetscViewer& viewer) const | |
876 | { | ||
877 | 8 | PetscErrorCode code = VecView(m_self,viewer); | |
878 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 8 times.
|
8 | CHKERRQ(code); |
879 | 8 | return code; | |
880 | } | ||
881 | |||
882 | /***********************************************************************************/ | ||
883 | /***********************************************************************************/ | ||
884 | |||
885 | // TODO: https://www.mcs.anl.gov/petsc/petsc-3.12/docs/manualpages/Viewer/PetscViewerFormat.html | ||
886 | |||
887 | 8 | inline PetscErrorCode PetscVector::saveInBinaryFormat(MPI_Comm comm, const std::string filename, const std::string folder) const | |
888 | { | ||
889 |
2/4✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 8 times.
✗ Branch 5 not taken.
|
8 | const std::string filenameBin = folder + filename + ".mb"; |
890 | PetscViewer binaryViewer; | ||
891 | PetscErrorCode code; | ||
892 |
2/6✓ Branch 2 taken 8 times.
✗ Branch 3 not taken.
✗ Branch 4 not taken.
✓ Branch 5 taken 8 times.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
|
8 | code = PetscViewerBinaryOpen(comm, filenameBin.c_str(), FILE_MODE_WRITE, &binaryViewer); CHKERRQ(code); |
893 |
2/6✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
8 | code = PetscViewerPushFormat(binaryViewer,PETSC_VIEWER_NATIVE); CHKERRQ(code); |
894 |
2/6✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
8 | code = this->view(binaryViewer); CHKERRQ(code); |
895 |
2/6✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
✗ Branch 6 not taken.
✗ Branch 7 not taken.
|
8 | code = PetscViewerDestroy(&binaryViewer); CHKERRQ(code); |
896 |
2/4✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
✓ Branch 4 taken 8 times.
|
8 | if ( FelisceParam::verbose() > 1 ) { |
897 | ✗ | const std::string msg = "Petsc Vector saved in " + filename + ".mb in folder " + folder + "\n"; | |
898 | ✗ | PetscPrintf(comm,"%s",msg.c_str()); | |
899 | } | ||
900 | 8 | return code; | |
901 | 8 | } | |
902 | |||
903 | /***********************************************************************************/ | ||
904 | /***********************************************************************************/ | ||
905 | |||
906 | ✗ | inline PetscErrorCode PetscVector::saveInMatlabFormat(MPI_Comm comm, const std::string filename, const std::string folder) const | |
907 | { | ||
908 | ✗ | const std::string filenameMatlab = folder + filename + ".m"; | |
909 | PetscViewer matlabViewer; | ||
910 | PetscErrorCode code; | ||
911 | ✗ | code = PetscViewerCreate(comm,&matlabViewer); | |
912 | ✗ | CHKERRQ(code); | |
913 | ✗ | code = PetscViewerASCIIOpen(comm,filenameMatlab.c_str(),&matlabViewer); | |
914 | ✗ | CHKERRQ(code); | |
915 | ✗ | code = PetscViewerPushFormat(matlabViewer,PETSC_VIEWER_ASCII_MATLAB); | |
916 | ✗ | CHKERRQ(code); | |
917 | ✗ | code = this->view(matlabViewer); | |
918 | ✗ | CHKERRQ(code); | |
919 | ✗ | code = PetscViewerDestroy(&matlabViewer); | |
920 | ✗ | CHKERRQ(code); | |
921 | ✗ | if ( FelisceParam::verbose() > 1 ) { | |
922 | ✗ | const std::string msg = "Petsc Vector saved in " + filename + ".m in folder " + folder + "\n"; | |
923 | ✗ | PetscPrintf(comm,"%s",msg.c_str()); | |
924 | } | ||
925 | ✗ | return code; | |
926 | } | ||
927 | |||
928 | //================ | ||
929 | // INPUT FUNCTIONS | ||
930 | //================ | ||
931 | |||
932 | ✗ | inline PetscErrorCode PetscVector::load(PetscViewer viewer) { | |
933 | ✗ | const PetscErrorCode code = VecLoad(m_self,viewer); | |
934 | ✗ | CHKERRQ(code); | |
935 | ✗ | return code; | |
936 | } | ||
937 | |||
938 | /***********************************************************************************/ | ||
939 | /***********************************************************************************/ | ||
940 | |||
941 | ✗ | inline PetscErrorCode PetscVector::loadFromBinaryFormat(MPI_Comm comm, const std::string filename, const std::string folder) | |
942 | { | ||
943 | ✗ | const std::string filenameBin = folder + filename + ".mb"; | |
944 | PetscViewer binaryViewer; | ||
945 | PetscErrorCode code; | ||
946 | ✗ | code = PetscViewerBinaryOpen(comm, filenameBin.c_str(), FILE_MODE_READ, &binaryViewer); CHKERRQ(code); | |
947 | ✗ | code = PetscViewerPushFormat(binaryViewer,PETSC_VIEWER_NATIVE); CHKERRQ(code); | |
948 | ✗ | code = this->load(binaryViewer); CHKERRQ(code); | |
949 | ✗ | code = PetscViewerDestroy(&binaryViewer); CHKERRQ(code); | |
950 | ✗ | if ( FelisceParam::verbose() > 1 ) { | |
951 | ✗ | const std::string msg = "Petsc Vector loaded from " + filename + ".mb in folder " + folder + "\n"; | |
952 | ✗ | PetscPrintf(comm,"%s",msg.c_str()); | |
953 | } | ||
954 | ✗ | return code; | |
955 | } | ||
956 | |||
957 | /***********************************************************************************/ | ||
958 | /***********************************************************************************/ | ||
959 | |||
960 | inline PetscErrorCode PetscVector::loadFromMatlabFormat(MPI_Comm comm, const std::string filename, const std::string folder) { | ||
961 | std::stringstream filenameBin; | ||
962 | filenameBin<< folder << filename << ".m"; | ||
963 | PetscViewer matlabViewer; | ||
964 | PetscErrorCode code; | ||
965 | // code = | ||
966 | PetscViewerMatlabOpen(comm, filenameBin.str().c_str(), FILE_MODE_READ, &matlabViewer); | ||
967 | // CHKERRQ(code); | ||
968 | code = PetscViewerPushFormat(matlabViewer,PETSC_VIEWER_ASCII_MATLAB); | ||
969 | CHKERRQ(code); | ||
970 | code = this->load(matlabViewer); | ||
971 | CHKERRQ(code); | ||
972 | code = PetscViewerDestroy(&matlabViewer); | ||
973 | CHKERRQ(code); | ||
974 | if ( FelisceParam::verbose() > 1 ) { | ||
975 | std::stringstream msg; | ||
976 | msg<<"Petsc Vector loaded from "<<filename<<".m in folder "<<folder<<std::endl; | ||
977 | PetscPrintf(comm,"%s",msg.str().c_str()); | ||
978 | } | ||
979 | return code; | ||
980 | } | ||
981 | |||
982 | /***********************************************************************************/ | ||
983 | /***********************************************************************************/ | ||
984 | |||
985 | /* Load .vct file */ | ||
986 | inline PetscErrorCode PetscVector::loadFromEnsight(const std::string fileName) | ||
987 | { | ||
988 | PetscErrorCode code; | ||
989 | |||
990 | felInt numDof; | ||
991 | code = this->getSize(&numDof); | ||
992 | PetscReal array[numDof]; | ||
993 | |||
994 | std::ifstream ensightFile ((fileName).c_str()); | ||
995 | std::string line; | ||
996 | |||
997 | getline(ensightFile,line); /* discard first line */ | ||
998 | PetscInt arrIndex = 0; | ||
999 | if (ensightFile.is_open()){ | ||
1000 | while (getline(ensightFile,line)){ | ||
1001 | if (line.length() == 72){ | ||
1002 | for (int i = 0; i < 6; i++){ | ||
1003 | array[arrIndex] = atof(line.substr(0, 12).c_str()); arrIndex++; | ||
1004 | line = line.substr(12, line.length()); | ||
1005 | } | ||
1006 | } else { | ||
1007 | for (int i = 0; i < 3; i++){ | ||
1008 | array[arrIndex] = atof(line.substr(0, 12).c_str()); arrIndex++; | ||
1009 | line = line.substr(12, line.length()); | ||
1010 | } | ||
1011 | } | ||
1012 | } | ||
1013 | } else { | ||
1014 | std::cout << "The file " + fileName + " cannot be opened \n"; | ||
1015 | } | ||
1016 | |||
1017 | PetscInt indexes[arrIndex + 1]; | ||
1018 | for (PetscInt i = 0; i < arrIndex + 1; i++) | ||
1019 | indexes[i] = i; | ||
1020 | |||
1021 | code = this->setValues(arrIndex + 1, indexes, array, INSERT_VALUES); | ||
1022 | CHKERRQ(code); | ||
1023 | |||
1024 | return code; | ||
1025 | } | ||
1026 | |||
1027 | //================= | ||
1028 | // TO BE DEPRECATED | ||
1029 | //================= | ||
1030 | |||
1031 | inline const Vec& PetscVector::toPetsc() const | ||
1032 | { | ||
1033 | return m_self; | ||
1034 | } | ||
1035 | |||
1036 | /***********************************************************************************/ | ||
1037 | /***********************************************************************************/ | ||
1038 | |||
1039 | 76977 | inline Vec& PetscVector::toPetsc() { | |
1040 | 76977 | return m_self; | |
1041 | } | ||
1042 | |||
1043 | /** | ||
1044 | * @brief Multiply (*) operator | ||
1045 | * @param scalar The value that multiplies the vector | ||
1046 | * @param b The vector to be multiplied | ||
1047 | * @return The dot product of the current vector and b | ||
1048 | */ | ||
1049 | 8 | inline PetscVector operator*( | |
1050 | const PetscInt scalar, | ||
1051 | const PetscVector& b | ||
1052 | ) | ||
1053 | { | ||
1054 | 8 | PetscVector a(b); | |
1055 |
1/2✓ Branch 1 taken 8 times.
✗ Branch 2 not taken.
|
8 | a.scale(PetscScalar(scalar)); |
1056 | 8 | return a; | |
1057 | } | ||
1058 | |||
1059 | /** | ||
1060 | * @brief Multiply (*) operator | ||
1061 | * @param b The vector to be multiplied | ||
1062 | * @param scalar The value that multiplies the vector | ||
1063 | * @return The dot product of the current vector and b | ||
1064 | */ | ||
1065 | 12 | inline PetscVector operator*( | |
1066 | const PetscScalar scalar, | ||
1067 | const PetscVector& b | ||
1068 | ) | ||
1069 | { | ||
1070 | 12 | PetscVector a(b); | |
1071 |
1/2✓ Branch 1 taken 12 times.
✗ Branch 2 not taken.
|
12 | a.scale(scalar); |
1072 | 12 | return a; | |
1073 | } | ||
1074 | |||
1075 | } | ||
1076 |