v0.15.5
Loading...
Searching...
No Matches
Schur.cpp
Go to the documentation of this file.
1/**
2 * @file Schur.cpp
3 * @brief Implementation of Schur Complement
4 * @date 2023-02-03
5 *
6 * @copyright Copyright (c) 2023
7 *
8 */
9
10namespace MoFEM {
11
12constexpr bool debug_schur = false;
13
16
17protected:
19
20 /**
21 * @brief Assemble Schur complement
22 *
23 */
25
26 Mat S,
27
28 const UId &uid_row, VectorInt &row_ind,
29
30 const UId &uid_col, VectorInt &col_ind,
31
32 MatrixDouble &m, InsertMode iora
33
34 );
35};
36
37// Assemble specialisations
39
40 Mat S,
41
42 const UId &uid_row, VectorInt &row_ind,
43
44 const UId &uid_col, VectorInt &col_ind,
45
46 MatrixDouble &m, InsertMode iora
47
48) {
49 return ::MatSetValues(
50
51 S, row_ind.size(), &*row_ind.begin(), col_ind.size(), &*col_ind.begin(),
52 &*m.data().begin(), iora
53
54 );
55}
56
57/**
58 * @brief Clear Schur complement internal data
59 *
60 */
62
64
65 MoFEMErrorCode doWork(int side, EntityType type,
67};
68
70
72 boost::shared_ptr<std::vector<unsigned char>> marker_ptr,
73 double diag_val);
74
75 MoFEMErrorCode doWork(int side, EntityType type,
77
78private:
79 boost::shared_ptr<std::vector<unsigned char>> markerPtr;
80 double diagVal;
81};
82
83/**
84 * @brief Assemble Schur complement (Implementation)
85 *
86 */
87template <typename OP_SCHUR_ASSEMBLE_BASE>
89
91
92 /**
93 * @brief Construct a new Op Schur Assemble End object
94 *
95 * @param fields_name list of fields
96 * @param field_ents list of entities on which schur complement is applied
97 * (can be empty)
98 * @param schur_ao map indices to schur matrix indices
99 * @param schur_mat schur matrix
100 * @param sym_schur true if schur is symmetric
101 * @param symm_op true if block diagonal is symmetric
102 */
104
105 std::vector<std::string> fields_name,
106 std::vector<boost::shared_ptr<Range>> field_ents,
107
108 SmartPetscObj<AO> schur_ao, SmartPetscObj<Mat> schur_mat,
109
110 bool sym_schur, bool symm_op);
111
112protected:
113 template <typename I>
115
116 int side, EntityType type, EntitiesFieldData::EntData &data
117
118 );
119
120 std::vector<std::string> fieldsName;
121 std::vector<boost::shared_ptr<Range>> fieldEnts;
125
129};
130
131struct SchurDSYSV; ///< SY symmetric
132struct SchurDGESV; ///< GE general (i.e., nonsymmetric, in some cases
133 ///< rectangular)
134
135/**
136 * @brief Assemble Schur complement
137 *
138 */
139template <typename I> struct OpSchurAssembleEnd;
140
141template <>
143 : public OpSchurAssembleEndImpl<OpSchurAssembleBaseImpl> {
144 using OpSchurAssembleEndImpl<OpSchurAssembleBaseImpl>::OpSchurAssembleEndImpl;
145 MoFEMErrorCode doWork(int side, EntityType type,
147};
148
149template <>
151 : public OpSchurAssembleEndImpl<OpSchurAssembleBaseImpl> {
152 using OpSchurAssembleEndImpl<OpSchurAssembleBaseImpl>::OpSchurAssembleEndImpl;
153 MoFEMErrorCode doWork(int side, EntityType type,
155};
156
157/**
158 * @brief Schur complement data storage
159 *
160 */
161struct SchurElemMats : public boost::enable_shared_from_this<SchurElemMats> {
162
163 SchurElemMats(const size_t idx, const UId uid_row, const UId uid_col);
164 virtual ~SchurElemMats() = default;
165
166 const size_t iDX;
169
170 inline auto &getMat() const { return locMats[iDX]; }
171 inline auto &getRowInd() const { return rowIndices[iDX]; }
172 inline auto &getColInd() const { return colIndices[iDX]; }
173
174 static MoFEMErrorCode MatSetValues(Mat M,
175 const EntitiesFieldData::EntData &row_data,
176 const EntitiesFieldData::EntData &col_data,
177 const MatrixDouble &mat, InsertMode iora);
178
179protected:
180 static MoFEMErrorCode
182 const EntitiesFieldData::EntData &col_data,
183 const MatrixDouble &mat, InsertMode iora);
184
185 friend struct OpSchurAssembleBegin;
187 template <typename OP_SCHUR_ASSEMBLE_BASE>
189
190 struct uid_mi_tag {};
191 struct col_mi_tag {};
192
193 using SchurElemStorage = multi_index_container<
194 const SchurElemMats *,
195 indexed_by<
196
197 ordered_unique<
198 tag<uid_mi_tag>,
199 composite_key<
201
202 member<SchurElemMats, const UId, &SchurElemMats::uidRow>,
203 member<SchurElemMats, const UId, &SchurElemMats::uidCol>
204
205 >>,
206
207 ordered_non_unique<tag<col_mi_tag>, member<SchurElemMats, const UId,
209
210 >>;
211
212 static boost::ptr_vector<MatrixDouble> locMats;
213 static boost::ptr_vector<VectorInt> rowIndices;
214 static boost::ptr_vector<VectorInt> colIndices;
215 static boost::ptr_vector<SchurElemMats> schurElemMats;
216 static size_t maxIndexCounter;
217
219};
220
222
223 virtual ~DiagBlockIndex() = default;
224
225 /**
226 * @brief block data indexes
227 *
228 */
229 struct Indexes {
230
236
237 inline UId getRowUId() const { return uid_row; }
238 inline UId getColUId() const { return uid_col; }
239 inline UId getFEUId() const { return uid_fe; }
240 inline int getRow() const { return row; }
241 inline int getCol() const { return col; }
242 inline int getNbRows() const { return nb_rows; }
243 inline int getNbCols() const { return nb_cols; }
244 inline int getLocRow() const { return loc_row; }
245 inline int getLocCol() const { return loc_col; }
246 inline int &getMatShift() const { return mat_shift; }
247
248 inline int rowShift() const {
249 return getRow() + getNbRows();
250 } // shift such that lower bound is included
251
252 inline int colShift() const {
253 return getCol() + getNbCols();
254 } // shift such that lower bound is included
255
256 private:
260 int row;
261 int col;
266 mutable int mat_shift;
267 };
268
269 using BlockIndex = multi_index_container<
270
271 Indexes,
272
273 indexed_by<
274
275 ordered_unique<
276
277 composite_key<Indexes,
278
279 const_mem_fun<Indexes, UId, &Indexes::getRowUId>,
280 const_mem_fun<Indexes, UId, &Indexes::getColUId>>
281
282 >,
283
284 ordered_non_unique<
285
286 const_mem_fun<Indexes, UId, &Indexes::getFEUId>
287
288 >,
289
290 ordered_non_unique<
291
292 const_mem_fun<Indexes, UId, &Indexes::getRowUId>
293
294 >,
295
296 ordered_non_unique<
297
298 const_mem_fun<Indexes, UId, &Indexes::getColUId>
299
300 >>>;
301
302 BlockIndex blockIndex; ///< blocks indexes storage
303};
304
306
309
310 boost::shared_ptr<std::vector<double>> dataBlocksPtr;
311 boost::shared_ptr<std::vector<double>> preconditionerBlocksPtr;
312 boost::shared_ptr<std::vector<double>> parentBlockStructurePtr;
313
315};
316
324
326 PetscLogEventRegister("schurMatSetVal", 0, &MOFEM_EVENT_schurMatSetValues);
327 PetscLogEventRegister("opSchurAsmEnd", 0, &MOFEM_EVENT_opSchurAssembleEnd);
328 PetscLogEventRegister("blockSetVal", 0, &MOFEM_EVENT_BlockStructureSetValues);
329 PetscLogEventRegister("blockMult", 0, &MOFEM_EVENT_BlockStructureMult);
330 PetscLogEventRegister("blockSolve", 0, &MOFEM_EVENT_BlockStructureSolve);
331 PetscLogEventRegister("schurZeroRandC", 0, &MOFEM_EVENT_zeroRowsAndCols);
332 PetscLogEventRegister("assembleSchurMat", 0, &MOFEM_EVENT_AssembleSchurMat);
333}
334
336boost::ptr_vector<MatrixDouble> SchurElemMats::locMats;
337boost::ptr_vector<VectorInt> SchurElemMats::rowIndices;
338boost::ptr_vector<VectorInt> SchurElemMats::colIndices;
339boost::ptr_vector<SchurElemMats> SchurElemMats::schurElemMats;
341
342SchurElemMats::SchurElemMats(const size_t idx, const UId uid_row,
343 const UId uid_col)
344 : iDX(idx), uidRow(uid_row), uidCol(uid_col) {}
345
348 const EntitiesFieldData::EntData &col_data,
349 const MatrixDouble &mat, InsertMode iora) {
351
352#ifndef NDEBUG
353 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_schurMatSetValues, 0, 0, 0, 0);
354#endif // NDEBUG
355
356 // get row indices, in case of store, get indices from storage
357 // storage keeps marked indices to manage boundary conditions
358 auto get_row_indices = [&]() -> const VectorInt & {
359 if (auto e_ptr = row_data.getFieldEntities()[0]) {
360 if (auto stored_data_ptr =
361 e_ptr->getSharedStoragePtr<EssentialBcStorage>()) {
362 return stored_data_ptr->entityIndices;
363 }
364 }
365 return row_data.getIndices();
366 };
367
368 const auto &row_ind = get_row_indices();
369 const auto &col_ind = col_data.getIndices();
370
371 const auto nb_rows = row_ind.size();
372 const auto nb_cols = col_ind.size();
373
374#ifndef NDEBUG
375 if (mat.size1() != nb_rows) {
376 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
377 "Wrong mat size %ld != %ld", mat.size1(), nb_rows);
378 }
379 if (mat.size2() != nb_cols) {
380 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
381 "Wrong mat size %ld != %ld", mat.size2(), nb_cols);
382 }
383#endif // NDEBUG
384
385 // get entity uid
386 auto get_uid = [](auto &data) {
387 if (data.getFieldEntities().size() == 1) {
388
389 return data.getFieldEntities()[0]->getLocalUniqueId();
390
391 } else {
392
393 // Is assumed that sum of entities ids gives unique id, that is not true,
394 // but corner case is improbable.
395
396 // @todo: debug version should have test
397
398 auto &uid0 = data.getFieldEntities()[0]->getLocalUniqueId();
399 auto field_id0 = FieldEntity::getFieldBitNumberFromUniqueId(uid0);
400 auto ent0 = FieldEntity::getHandleFromUniqueId(uid0);
401 auto type0 = type_from_handle(ent0);
402 auto id = id_from_handle(ent0);
403
404 for (auto i = 1; i < data.getFieldEntities().size(); ++i) {
405
406 // get entity id from ent
407 id += id_from_handle(
408
409 // get entity handle from unique uid
411 data.getFieldEntities()[i]->getLocalUniqueId())
412
413 );
414 }
415
417 field_id0,
418
419 ent_form_type_and_id(type0, id)
420
421 );
422 }
423 };
424
425 auto uid_row = get_uid(row_data);
426 auto uid_col = get_uid(col_data);
427
428 auto it =
429 SchurElemMats::schurL2Storage.template get<SchurElemMats::uid_mi_tag>()
430 .find(boost::make_tuple(uid_row, uid_col));
431
432 if (it ==
433 SchurElemMats::schurL2Storage.template get<SchurElemMats::uid_mi_tag>()
434 .end()) {
435
436 // get size of arrays of matrices
437 const auto size = SchurElemMats::locMats.size();
438
439 // expand memory allocation
440 if (SchurElemMats::maxIndexCounter == size) {
441 SchurElemMats::locMats.push_back(new MatrixDouble());
442 SchurElemMats::rowIndices.push_back(new VectorInt());
443 SchurElemMats::colIndices.push_back(new VectorInt());
445 new SchurElemMats(SchurElemMats::maxIndexCounter, uid_row, uid_col));
446 } else {
448 uid_row;
450 uid_col;
451 }
452
453 // add matrix to storage
454 auto p = SchurElemMats::schurL2Storage.emplace(
456#ifndef NDEBUG
457 if (!p.second) {
458 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Failed to insert");
459 }
460#endif // NDEBUG
461
462 auto asmb = [&](auto &sm) {
463 sm.resize(nb_rows, nb_cols, false);
464 noalias(sm) = mat;
465 };
466
467 asmb((*p.first)->getMat());
468
469 auto add_indices = [](auto &storage, auto &ind) {
470 storage.resize(ind.size(), false);
471 noalias(storage) = ind;
472 };
473
474 add_indices((*p.first)->getRowInd(), row_ind);
475 add_indices((*p.first)->getColInd(), col_ind);
476
477 } else {
478 // entry (submatrix) already exists
479
480 auto asmb = [&](auto &sm) {
482
483#ifndef NDEBUG
484 if (sm.size1() != nb_rows) {
485 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
486 "Wrong mat or storage size %ld != %ld", sm.size1(), nb_rows);
487 }
488 if (sm.size2() != nb_cols) {
489 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
490 "Wrong mat or storage size %ld != %ld", sm.size2(), nb_cols);
491 }
492#endif // NDEBUG
493
494 switch (iora) {
495 case ADD_VALUES:
496 sm += mat;
497 break;
498 case INSERT_VALUES:
499 noalias(sm) = mat;
500 break;
501 default:
502 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED,
503 "Assembly type not implemented");
504 }
506 };
507
508 CHKERR asmb((*it)->getMat());
509
510 // no need to set indices
511 }
512
513#ifndef NDEBUG
514 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_schurMatSetValues, 0, 0, 0, 0);
515#endif // NDEBUG
516
518}
519
522
526#ifndef NDEBUG
527 if constexpr (debug_schur)
528 MOFEM_LOG("SELF", Sev::noisy) << "Schur assemble begin";
529#endif
532
534}
535
537 boost::shared_ptr<std::vector<unsigned char>> marker_ptr, double diag_val)
538 : OpSchurAssembleBaseImpl(NOSPACE, OPSPACE), markerPtr(marker_ptr),
539 diagVal(diag_val) {}
540
542OpSchurZeroRowsAndCols::doWork(int side, EntityType type,
545
546 auto fe_ptr = getFEMethod();
547 auto prb_ptr = fe_ptr->problemPtr;
548 auto &dof = prb_ptr->getNumeredRowDofsPtr()->get<PetscGlobalIdx_mi_tag>();
549
550 auto find_loc_dof_index = [&](auto glob) {
551 auto it = dof.find(glob);
552 if (it != dof.end()) {
553 return (*it)->getPetscLocalDofIdx();
554 }
555 MOFEM_LOG("SELF", Sev::error) << "Global dof " << glob;
558 "Global dof index not found in local numered dofs storage");
559 return -1;
560 };
561
562 for (auto &s : SchurElemMats::schurL2Storage) {
563
564 auto zero_ind = [&](auto &ind) {
565 for (auto i = 0; i < ind.size(); ++i) {
566 if (ind[i] >= 0) {
567 if ((*markerPtr)[find_loc_dof_index(ind[i])]) {
568 ind[i] = -1;
569 }
570 }
571 }
572 };
573
574 zero_ind(s->getRowInd());
575 zero_ind(s->getColInd());
576
577 }
578
580};
581
582template <typename OP_SCHUR_ASSEMBLE_BASE>
584 std::vector<std::string> fields_name,
585 std::vector<boost::shared_ptr<Range>> field_ents,
586 SmartPetscObj<AO> schur_ao, SmartPetscObj<Mat> schur_mat, bool sym_schur,
587 bool symm_op)
588 : OP(NOSPACE, OP::OPSPACE, symm_op), fieldsName(fields_name),
589 fieldEnts(field_ents), schurAO(schur_ao), schurMat(schur_mat),
590 symSchur(sym_schur) {}
591
592template <typename OP_SCHUR_ASSEMBLE_BASE>
593template <typename I>
595 int side, EntityType type, EntitiesFieldData::EntData &data) {
596
598
599#ifndef NDEBUG
600 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_opSchurAssembleEnd, 0, 0, 0, 0);
601#endif
602
603#ifndef NDEBUG
604 if constexpr (debug_schur)
605 MOFEM_LOG("SELF", Sev::noisy) << "Schur assemble begin -> end";
606#endif
607
608 auto get_a00_uids = [&]() {
609 auto get_field_bit = [&](auto &name) {
610 return OP::getPtrFE()->mField.get_field_bit_number(name);
611 };
612
613 std::vector<std::pair<UId, UId>> a00_uids;
614 a00_uids.reserve(fieldsName.size());
615 for (auto ss = 0; ss != fieldsName.size(); ++ss) {
616 auto field_bit = get_field_bit(fieldsName[ss]);
617 auto row_ents = fieldEnts[ss];
618 if (row_ents) {
619 for (auto p = row_ents->pair_begin(); p != row_ents->pair_end(); ++p) {
620 auto lo_uid =
621 FieldEntity::getLoLocalEntityBitNumber(field_bit, p->first);
622 auto hi_uid =
623 FieldEntity::getHiLocalEntityBitNumber(field_bit, p->second);
624 a00_uids.push_back(std::make_pair(lo_uid, hi_uid));
625 }
626 } else {
628 field_bit, get_id_for_min_type<MBVERTEX>());
630 field_bit, get_id_for_max_type<MBENTITYSET>());
631 a00_uids.push_back(std::make_pair(lo_uid, hi_uid));
632 }
633 }
634 return a00_uids;
635 };
636
637#ifndef NDEBUG
638 auto get_field_name = [&](auto uid) {
639 return OP::getPtrFE()->mField.get_field_name(field_bit_from_bit_number(
641 };
642 auto list_storage = [&](auto &storage) {
644 int i = 0;
645 for (auto &p : storage) {
646 MOFEM_LOG("SELF", Sev::noisy)
647 << "List schur storage: " << i << " " << p->iDX << ": "
648 << get_field_name(p->uidRow) << " " << get_field_name(p->uidCol)
649 << " : " << p->getMat().size1() << " " << p->getMat().size2();
650 ++i;
651 }
653 };
654#endif // NDEBUG
655
656 auto assemble_dense_blocks = [&]() {
657 using matrix_range = ublas::matrix_range<MatrixDouble>;
658 using range = ublas::range;
660 auto &storage = SchurElemMats::schurL2Storage;
661
662 auto assemble_schur = [this](auto &m, auto &uid_row, auto &uid_col,
663 auto *row_ind_ptr, auto *col_ind_ptr) {
665
666 if (schurMat) {
667
668 if (auto ierr = this->assembleSchurMat(
669
670 schurMat, uid_row, *row_ind_ptr, uid_col, *col_ind_ptr, m,
671 ADD_VALUES
672
673 )) {
674#ifndef NDEBUG
675 auto field_ents = OP::getPtrFE()->mField.get_field_ents();
676 auto row_ent_it = field_ents->find(uid_row);
677 auto col_ent_it = field_ents->find(uid_col);
678 MOFEM_LOG_CHANNEL("SELF");
679 if (row_ent_it != field_ents->end())
680 MOFEM_LOG("SELF", Sev::error)
681 << "Assemble row entity: " << (*row_ent_it)->getName() << " "
682 << (*col_ent_it)->getEntTypeName() << " side "
683 << (*row_ent_it)->getSideNumber();
684 if (col_ent_it != field_ents->end())
685 MOFEM_LOG("SELF", Sev::error)
686 << "Assemble col entity: " << (*col_ent_it)->getName() << " "
687 << (*col_ent_it)->getEntTypeName() << " side "
688 << (*col_ent_it)->getSideNumber();
689#endif // NDEBUG
690 CHK_THROW_MESSAGE(ierr, "MatSetValues");
691 }
692 }
693
695 };
696
697 auto a00_uids = get_a00_uids();
698
699 auto get_block_indexing = [&](auto &a00_uids) {
700 // iterate over a00 uids and find blocks
701 std::vector<const SchurElemMats *> block_list;
702 block_list.reserve(storage.size());
703 for (auto &rp_uid : a00_uids) {
704 auto [rlo_uid, rhi_uid] = rp_uid;
705 for (auto &cp_uid : a00_uids) {
706 auto [clo_uid, chi_uid] = cp_uid;
707
708 auto it =
709 storage.template get<SchurElemMats::uid_mi_tag>().lower_bound(
710 boost::make_tuple(rlo_uid, clo_uid));
711 auto hi_it =
712 storage.template get<SchurElemMats::uid_mi_tag>().upper_bound(
713 boost::make_tuple(rhi_uid, chi_uid));
714
715 for (; it != hi_it; ++it) {
716 if ((*it)->uidRow >= rlo_uid && (*it)->uidRow < rhi_uid &&
717 (*it)->uidCol >= clo_uid && (*it)->uidCol < chi_uid) {
718 block_list.push_back(*it);
719 }
720 }
721 }
722 }
723
724 // create block indexes map for blockMat
725 std::map<UId, std::pair<size_t, const VectorInt *>>
726 block_indexing; // uid block map
727 for (auto d : block_list) {
728 if (block_indexing.find(d->uidRow) == block_indexing.end()) {
729 block_indexing[d->uidRow] =
730 std::make_pair(d->getRowInd().size(), &(d->getRowInd()));
731 }
732 }
733
734 // set indexes to block
735 int mat_block_size = 0; // size of block matrix
736 for (auto &p_uid : a00_uids) {
737 auto [lo_uid, hi_uid] = p_uid;
738 auto lo = block_indexing.lower_bound(lo_uid);
739 auto up = block_indexing.upper_bound(hi_uid);
740 for (; lo != up; ++lo) {
741 lo->second.first = mat_block_size;
742 mat_block_size += lo->second.second->size();
743 }
744 }
745
746 return std::make_tuple(block_list, block_indexing, mat_block_size);
747 };
748
749 auto get_schur_block_list = [&](auto &block_indexing) {
750 std::vector<const SchurElemMats *> block_list;
751 block_list.reserve(storage.size());
752 for (auto &s : storage) {
753 if (block_indexing.find(s->uidRow) == block_indexing.end() &&
754 block_indexing.find(s->uidCol) == block_indexing.end()) {
755 block_list.push_back(s);
756 }
757 }
758 return block_list;
759 };
760
761 auto bi = get_block_indexing(a00_uids);
762 auto &[block_list, block_indexing, block_mat_size] = bi;
763
764 if (block_mat_size == 0) {
765
766 if (schurAO) {
767 for (auto &s : storage) {
768 CHKERR AOApplicationToPetsc(schurAO, s->getRowInd().size(),
769 &*s->getRowInd().begin());
770 CHKERR AOApplicationToPetsc(schurAO, s->getColInd().size(),
771 &*s->getColInd().begin());
772 }
773 }
774
775 for (auto &s : storage) {
776 auto &m = s->getMat();
777 CHKERR assemble_schur(m, s->uidRow, s->uidCol, &(s->getRowInd()),
778 &(s->getColInd()));
779 }
781 }
782
783 blockMat.resize(block_mat_size, block_mat_size, false);
784 blockMat.clear();
785
786 auto get_range = [](auto &bi) {
787 return range(bi.first, bi.first + bi.second->size());
788 };
789
790 for (auto &s : block_list) {
791 auto &m = s->getMat();
792#ifndef NDEBUG
793 if (block_indexing.find(s->uidRow) == block_indexing.end())
794 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong rlo_uid");
795 if (block_indexing.find(s->uidCol) == block_indexing.end())
796 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "Wrong clo_uid");
797#endif // NDEBUG
798
799 auto &rbi = block_indexing.at(s->uidRow);
800 auto &cbi = block_indexing.at(s->uidCol);
801
802 auto sub_mat = matrix_range(blockMat, get_range(rbi), get_range(cbi));
803 sub_mat = m;
804 }
805
806 auto get_zeroed_indices = [&](auto extractor_uid, auto extractor_ind) {
807 auto &[block_list, block_indexing, block_mat_size] = bi;
808 std::vector<int> zeroed_indices;
809 zeroed_indices.reserve(block_mat_size);
810 for (auto &s : block_list) {
811 auto &bi = block_indexing.at(extractor_uid(s));
812 auto &ind = extractor_ind(s);
813 for (auto it = ind.begin(); it != ind.end(); ++it) {
814 if (*it < 0) {
815 auto idx = bi.first + std::distance(ind.begin(), it);
816 zeroed_indices.push_back(idx);
817 }
818 }
819 }
820 std::sort(zeroed_indices.begin(), zeroed_indices.end());
821 auto it = std::unique(zeroed_indices.begin(), zeroed_indices.end());
822 zeroed_indices.resize(std::distance(zeroed_indices.begin(), it));
823 return zeroed_indices;
824 };
825 auto zero_rows = get_zeroed_indices(
826 [](auto &s) { return s->uidRow; },
827 [](auto &s) -> VectorInt & { return s->getRowInd(); });
828 auto zero_cols = get_zeroed_indices(
829 [](auto &s) { return s->uidCol; },
830 [](auto &s) -> VectorInt & { return s->getColInd(); });
831
832 for (auto i : zero_rows) {
833 for (auto j = 0; j != blockMat.size2(); ++j) {
834 blockMat(i, j) = 0;
835 }
836 }
837 for (auto j : zero_cols) {
838 for (auto i = 0; i != blockMat.size1(); ++i) {
839 blockMat(i, j) = 0;
840 }
841 }
842 for (auto i : zero_rows) {
843 blockMat(i, i) = 1;
844 }
845
846 CHKERR I::invertMat(blockMat, invMat);
847
848 // clear storage and block list from a00 blocks, no more needed
849 for (auto &s : block_list) {
850 auto it = storage.template get<SchurElemMats::uid_mi_tag>().find(
851 boost::make_tuple(s->uidRow, s->uidCol));
852 storage.template get<SchurElemMats::uid_mi_tag>().erase(it);
853 }
854 block_list.clear();
855
856 if (schurAO) {
857 for (auto &s : storage) {
858 if (block_indexing.find(s->uidRow) == block_indexing.end()) {
859 CHKERR AOApplicationToPetsc(schurAO, s->getRowInd().size(),
860 &*s->getRowInd().begin());
861 }
862 if (block_indexing.find(s->uidCol) == block_indexing.end()) {
863 CHKERR AOApplicationToPetsc(schurAO, s->getColInd().size(),
864 &*s->getColInd().begin());
865 }
866 }
867 }
868
869 auto schur_block_list = get_schur_block_list(block_indexing);
870 for (auto &s : schur_block_list) {
871 auto &m = s->getMat();
872 CHKERR assemble_schur(m, s->uidRow, s->uidCol, &(s->getRowInd()),
873 &(s->getColInd()));
874 }
875 for (auto &s : schur_block_list) {
876 auto it = storage.template get<SchurElemMats::uid_mi_tag>().find(
877 boost::make_tuple(s->uidRow, s->uidCol));
878 storage.template get<SchurElemMats::uid_mi_tag>().erase(it);
879 }
880 schur_block_list.clear();
881
882 for (
883
884 auto rp_uid_it = a00_uids.begin(); rp_uid_it != a00_uids.end();
885 ++rp_uid_it
886
887 ) {
888 auto [rlo_uid, rhi_uid] = *rp_uid_it;
889 for (auto rm = block_indexing.lower_bound(rlo_uid);
890 rm != block_indexing.upper_bound(rhi_uid); ++rm) {
891 auto &rbi = rm->second;
892
893 auto a_lo_tmp =
894 storage.template get<SchurElemMats::col_mi_tag>().lower_bound(
895 rm->first);
896 auto a_hi =
897 storage.template get<SchurElemMats::col_mi_tag>().upper_bound(
898 rm->first);
899
900 for (
901
902 auto cp_uid_it = (symSchur) ? rp_uid_it : a00_uids.begin();
903 cp_uid_it != a00_uids.end(); ++cp_uid_it
904
905 ) {
906 auto [clo_uid, chi_uid] = *cp_uid_it;
907 for (auto cm = block_indexing.lower_bound(clo_uid);
908 cm != block_indexing.upper_bound(chi_uid); ++cm) {
909 auto &cbi = cm->second;
910
911 auto c_lo_tmp =
912 storage.template get<SchurElemMats::uid_mi_tag>().lower_bound(
913 boost::make_tuple(cm->first, 0));
914 auto c_hi =
915 storage.template get<SchurElemMats::uid_mi_tag>().upper_bound(
916 boost::make_tuple(cm->first, FieldEntity::getHiBitNumberUId(
917 BITFIELDID_SIZE - 1)));
918
919 auto sub_inv_mat =
920 matrix_range(invMat, get_range(rbi), get_range(cbi));
921 bM.resize(sub_inv_mat.size1(), sub_inv_mat.size2());
922 noalias(bM) = sub_inv_mat;
923
924 for (auto a_lo = a_lo_tmp; a_lo != a_hi; ++a_lo) {
925#ifndef NDEBUG
926 if (block_indexing.find((*a_lo)->uidRow) != block_indexing.end())
927 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
928 "Wrong a_lo->uidRow");
929#endif
930
931 auto &a = (*a_lo)->getMat();
932 abM.resize(a.size1(), bM.size2(), false);
933 abM.clear();
934 cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans,
935
936 a.size1(), bM.size2(), a.size2(), 1.,
937
938 &*a.data().begin(), a.size2(),
939
940 &*bM.data().begin(), bM.size2(), 0.,
941
942 &*abM.data().begin(), abM.size2());
943
944 for (auto c_lo = c_lo_tmp; c_lo != c_hi; ++c_lo) {
945#ifndef NDEBUG
946 if (block_indexing.find((*c_lo)->uidCol) !=
947 block_indexing.end())
948 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
949 "Wrong a_lo->uidRow");
950#endif
951
952 auto &c = (*c_lo)->getMat();
953
954 abcM.resize(abM.size1(), c.size2(), false);
955 abcM.clear();
956
957 cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasNoTrans,
958 abM.size1(), c.size2(), abM.size2(), -1.,
959
960 &*abM.data().begin(), abM.size2(),
961
962 &*c.data().begin(), c.size2(), 0.,
963
964 &*abcM.data().begin(), abcM.size2());
965
966 CHKERR assemble_schur(abcM, (*a_lo)->uidRow, (*c_lo)->uidCol,
967 &((*a_lo)->getRowInd()),
968 &((*c_lo)->getColInd()));
969
970 if (symSchur && rp_uid_it != cp_uid_it) {
971 abcM = trans(abcM);
972 CHKERR assemble_schur(abcM, (*c_lo)->uidCol, (*a_lo)->uidRow,
973 &((*c_lo)->getColInd()),
974 &((*a_lo)->getRowInd()));
975 }
976 }
977 }
978 }
979 }
980 }
981 }
982
984 };
985
986 // Assemble Schur complements
987 CHKERR assemble_dense_blocks();
988
989#ifndef NDEBUG
990 if constexpr (debug_schur)
991 MOFEM_LOG("SELF", Sev::noisy) << "Schur assemble done";
992#endif
993
994#ifndef NDEBUG
995 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_opSchurAssembleEnd, 0, 0, 0, 0);
996#endif
997
999}
1000
1002
1003 static auto invertMat(MatrixDouble &m, MatrixDouble &inv) {
1005
1006 const auto nb = m.size1();
1007#ifndef NDEBUG
1008 if (nb != m.size2()) {
1009 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1010 "It should be square matrix %ld != %ld", nb, m.size2());
1011 }
1012#endif
1013
1014 inv.resize(nb, nb, false);
1015 inv.swap(m);
1016 m.resize(2 * nb, 2 * nb, false);
1017
1018 VectorInt ipiv(nb);
1019 int info;
1020
1021 // FIXME: add lapack_dsytrf and lapack_dsytri
1022
1023 info = lapack_dgetrf(nb, nb, &*inv.data().begin(), nb, &*ipiv.begin());
1024 if (info)
1025 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1026 "lapack error info = %d", info);
1027 info = lapack_dgetri(nb, &*inv.data().begin(), nb, &*ipiv.begin(),
1028 &*m.data().begin(), m.data().size());
1029 if (info)
1030 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1031 "lapack error info = %d", info);
1032
1034 }
1035};
1036
1038
1039 static auto invertMat(MatrixDouble &m, MatrixDouble &inv) {
1041
1042 const auto nb = m.size1();
1043#ifndef NDEBUG
1044 if (nb != m.size2()) {
1045 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1046 "It should be square matrix %ld != %ld", nb, m.size2());
1047 }
1048#endif
1049
1050 inv.resize(nb, nb, false);
1051 inv.swap(m);
1052
1053 VectorInt ipiv(nb);
1054 int info;
1055
1056 info = lapack_dgetrf(nb, nb, &*inv.data().begin(), nb, &*ipiv.begin());
1057 if (info)
1058 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1059 "lapack error info = %d", info);
1060 info = lapack_dgetri(nb, &*inv.data().begin(), nb, &*ipiv.begin(),
1061 &*m.data().begin(), m.data().size());
1062 if (info)
1063 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1064 "lapack error info = %d", info);
1065
1067 }
1068};
1069
1073 return doWorkImpl<SchurDSYSV>(side, type, data);
1074}
1075
1079 return doWorkImpl<SchurDGESV>(side, type, data);
1080}
1081
1082boost::shared_ptr<BlockStructure> createBlockMatStructure(
1083
1084 DM dm, //< dm
1085 SchurFEOpsFEandFields schur_fe_op_vec //< block elements
1086
1087) {
1088
1089 auto cmp_uid_lo = [](const boost::weak_ptr<FieldEntity> &a, const UId &b) {
1090 if (auto a_ptr = a.lock()) {
1091 if (a_ptr->getLocalUniqueId() < b)
1092 return true;
1093 else
1094 return false;
1095 } else {
1096 return false;
1097 }
1098 };
1099
1100 auto cmp_uid_hi = [](const UId &b, const boost::weak_ptr<FieldEntity> &a) {
1101 if (auto a_ptr = a.lock()) {
1102 if (b < a_ptr->getLocalUniqueId())
1103 return true;
1104 else
1105 return false;
1106 } else {
1107 return true;
1108 }
1109 };
1110
1111 // get uids for fields
1112 auto get_uid_pair = [](const auto &field_id) {
1114 field_id, get_id_for_min_type<MBVERTEX>());
1116 field_id, get_id_for_max_type<MBENTITYSET>());
1117 return std::make_pair(lo_uid, hi_uid);
1118 };
1119
1120 // get uids pair
1121 auto get_it_pair = [cmp_uid_lo, cmp_uid_hi](auto &&field_ents, auto &&p_uid) {
1122 auto lo = std::lower_bound(field_ents.begin(), field_ents.end(),
1123 p_uid.first, cmp_uid_lo);
1124 auto hi = std::upper_bound(field_ents.begin(), field_ents.end(),
1125 p_uid.second, cmp_uid_hi);
1126 return std::make_pair(lo, hi);
1127 };
1128
1129 // extract DOFs for rows/columns. DOFs are associated with fields entities
1130 // for given problem.
1131 auto row_extractor = [](auto &e) { return e->entityCacheRowDofs; };
1132 auto col_extractor = [](auto &e) { return e->entityCacheColDofs; };
1133
1134 auto extract_data = [](auto &&its, auto extractor) {
1135 std::vector<std::tuple<UId, int, int, int>> data;
1136 data.reserve(std::distance(its.first, its.second));
1137 // iterate field dofs
1138 for (; its.first != its.second; ++its.first) {
1139 if (auto e = its.first->lock()) {
1140 if (auto cache = extractor(e).lock()) {
1141 int nb_dofs = 0;
1142 for (auto dof = cache->loHi[0]; dof != cache->loHi[1]; ++dof) {
1143 if ((*dof)->getPetscGlobalDofIdx() != -1)
1144 ++nb_dofs;
1145 }
1146 auto uid = e->getLocalUniqueId();
1147 if (nb_dofs) {
1148
1149 auto glob = (*cache->loHi[0])->getPetscGlobalDofIdx();
1150 auto loc = (*cache->loHi[0])->getPetscLocalDofIdx();
1151 while (glob == -1 && cache->loHi[0] != cache->loHi[1]) {
1152 ++cache->loHi[0];
1153 glob = (*cache->loHi[0])->getPetscGlobalDofIdx();
1154 loc = (*cache->loHi[0])->getPetscLocalDofIdx();
1155 }
1156 data.emplace_back(uid, glob, nb_dofs, loc);
1157
1158#ifndef NDEBUG
1159
1160 for (auto lo = cache->loHi[0]; lo != cache->loHi[1]; ++lo) {
1161 auto glob = (*lo)->getPetscGlobalDofIdx();
1162 if (glob == -1) {
1164 "Wrong global index");
1165 }
1166 }
1167
1168#endif
1169 } else {
1170 data.emplace_back(uid, -1, 0, -1);
1171 }
1172 }
1173 }
1174 }
1175 return data;
1176 };
1177
1178 auto data_ptr = boost::make_shared<BlockStructure>();
1179 auto m_field_ptr = getInterfacePtr(dm);
1180
1181 // create element to extract data
1182 auto fe_method = boost::shared_ptr<MoFEM::FEMethod>(new MoFEM::FEMethod());
1183
1184 for (auto &d : schur_fe_op_vec) {
1185
1186 // extract bit numbers for row and column
1187 auto get_bit_numbers = [&d](auto op) {
1188 std::vector<FieldBitNumber> bit_numbers(d.second.size());
1189 std::transform(d.second.begin(), d.second.end(), bit_numbers.begin(), op);
1190 return bit_numbers;
1191 };
1192
1193 // extract bit numbers for row
1194 auto row_bit_numbers = get_bit_numbers(
1195 [&](auto &p) { return m_field_ptr->get_field_bit_number(p.first); });
1196 // extract bit numbers for row
1197 auto col_bit_numbers = get_bit_numbers(
1198 [&](auto &p) { return m_field_ptr->get_field_bit_number(p.second); });
1199
1200 fe_method->preProcessHook = []() { return 0; };
1201 fe_method->postProcessHook = []() { return 0; };
1202 fe_method->operatorHook = [&]() {
1204
1205 auto fe_uid = fe_method->numeredEntFiniteElementPtr->getLocalUniqueId();
1206
1207 for (auto f = 0; f != row_bit_numbers.size(); ++f) {
1208
1209 auto row_data =
1210 extract_data(get_it_pair(fe_method->getRowFieldEnts(),
1211 get_uid_pair(row_bit_numbers[f])),
1212 row_extractor);
1213 auto col_data =
1214 extract_data(get_it_pair(fe_method->getColFieldEnts(),
1215 get_uid_pair(col_bit_numbers[f])),
1216 col_extractor);
1217
1218 for (auto &r : row_data) {
1219 auto [r_uid, r_glob, r_nb_dofs, r_loc] = r;
1220 for (auto &c : col_data) {
1221 auto [c_uid, c_glob, c_nb_dofs, c_loc] = c;
1222 data_ptr->blockIndex.insert(BlockStructure::Indexes{
1223 r_uid, c_uid, fe_uid, r_glob, c_glob, r_nb_dofs, c_nb_dofs,
1224 r_loc, c_loc, -1});
1225 }
1226 }
1227 }
1228
1230 };
1231
1232 CHKERR DMoFEMLoopFiniteElementsUpAndLowRank(dm, d.first, fe_method, 0,
1233 m_field_ptr->get_comm_size());
1234 };
1235
1236 // order by column (that is for matrix multiplication)
1237 auto mem_size = 0;
1238 for (auto &v : data_ptr->blockIndex.get<0>()) {
1239 v.getMatShift() = mem_size;
1240 mem_size += v.getNbCols() * v.getNbRows();
1241 }
1242
1243 std::vector<double> tmp;
1244 if (mem_size > tmp.max_size())
1246
1247 data_ptr->dataBlocksPtr =
1248 boost::make_shared<std::vector<double>>(mem_size, 0);
1249
1250 auto ghost_x = createDMVector(dm);
1251 auto ghost_y = createDMVector(dm);
1252 CHKERR VecSetDM(ghost_x, PETSC_NULLPTR);
1253 CHKERR VecSetDM(ghost_y, PETSC_NULLPTR);
1254
1255 data_ptr->ghostX = ghost_x;
1256 data_ptr->ghostY = ghost_y;
1257
1258 return data_ptr;
1259}
1260
1262 Mat mat, Vec x, Vec y, InsertMode iora,
1263 boost::function<
1264 int(DiagBlockIndex::BlockIndex::nth_index<0>::type::iterator)>
1265 shift_extractor,
1266 boost::shared_ptr<std::vector<double>> data_blocks_ptr,
1267 bool multiply_by_preconditioner, enum CBLAS_TRANSPOSE trans_A);
1268
1269static MoFEMErrorCode solve_schur_block_shell(Mat mat, Vec y, Vec x,
1270 InsertMode iora,
1271 enum CBLAS_TRANSPOSE trans_A);
1272
1273static PetscErrorCode mult(Mat mat, Vec x, Vec y) {
1274 BlockStructure *ctx;
1275 CHKERR MatShellGetContext(mat, (void **)&ctx);
1277 mat, x, y, INSERT_VALUES,
1278
1279 [](DiagBlockIndex::BlockIndex::nth_index<0>::type::iterator it) {
1280 return it->getMatShift();
1281 },
1282
1283 ctx->dataBlocksPtr, true, CblasNoTrans);
1284}
1285static PetscErrorCode mult_add(Mat mat, Vec x, Vec y) {
1286 BlockStructure *ctx;
1287 CHKERR MatShellGetContext(mat, (void **)&ctx);
1289 mat, x, y, ADD_VALUES,
1290
1291 [](DiagBlockIndex::BlockIndex::nth_index<0>::type::iterator it) {
1292 return it->getMatShift();
1293 },
1294
1295 ctx->dataBlocksPtr, true, CblasNoTrans);
1296}
1297static PetscErrorCode mult_transpose(Mat mat, Vec x, Vec y) {
1298 BlockStructure *ctx;
1299 CHKERR MatShellGetContext(mat, (void **)&ctx);
1301 mat, x, y, INSERT_VALUES,
1302
1303 [](DiagBlockIndex::BlockIndex::nth_index<0>::type::iterator it) {
1304 return it->getMatShift();
1305 },
1306
1307 ctx->dataBlocksPtr, true, CblasTrans);
1308}
1309static PetscErrorCode mult_transpose_add(Mat mat, Vec x, Vec y) {
1310 BlockStructure *ctx;
1311 CHKERR MatShellGetContext(mat, (void **)&ctx);
1313 mat, x, y, ADD_VALUES,
1314
1315 [](DiagBlockIndex::BlockIndex::nth_index<0>::type::iterator it) {
1316 return it->getMatShift();
1317 },
1318
1319 ctx->dataBlocksPtr, true, CblasTrans);
1320}
1321
1322static PetscErrorCode solve(Mat mat, Vec x, Vec y) {
1323 return solve_schur_block_shell(mat, x, y, INSERT_VALUES, CblasNoTrans);
1324}
1325static PetscErrorCode solve_add(Mat mat, Vec x, Vec y) {
1326 return solve_schur_block_shell(mat, x, y, ADD_VALUES, CblasNoTrans);
1327}
1328static PetscErrorCode solve_transpose(Mat mat, Vec x, Vec y) {
1329 return solve_schur_block_shell(mat, x, y, INSERT_VALUES, CblasTrans);
1330}
1331static PetscErrorCode solve_transpose_add(Mat mat, Vec x, Vec y) {
1332 return solve_schur_block_shell(mat, x, y, ADD_VALUES, CblasTrans);
1333}
1334
1335static PetscErrorCode zero_rows_columns(Mat A, PetscInt N,
1336 const PetscInt rows[], PetscScalar diag,
1337 Vec x, Vec b) {
1338
1340 BlockStructure *ctx;
1341 CHKERR MatShellGetContext(A, (void **)&ctx);
1342
1343 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_zeroRowsAndCols, 0, 0, 0, 0);
1344
1345 using BlockIndexView = multi_index_container<
1346
1348
1349 indexed_by<
1350
1351 ordered_non_unique<
1352
1353 const_mem_fun<DiagBlockIndex::Indexes, int,
1355
1356 >,
1357
1358 ordered_non_unique<
1359
1360 const_mem_fun<DiagBlockIndex::Indexes, int,
1362
1363 >
1364
1365 >>;
1366
1367 BlockIndexView view;
1368 auto hint = view.get<0>().end();
1369 for (auto &v : ctx->blockIndex) {
1370 hint = view.insert(hint, &v);
1371 }
1372
1373 const int *ptr = &rows[0];
1374 int is_nb_rows = N;
1375 SmartPetscObj<IS> is_local;
1376
1377 MPI_Comm comm;
1378 CHKERR PetscObjectGetComm((PetscObject)A, &comm);
1379 int size;
1380 MPI_Comm_size(comm, &size);
1381 if (size > 1) {
1382 auto is = createISGeneral(comm, N, rows, PETSC_USE_POINTER);
1383 is_local = isAllGather(is);
1384 }
1385 if (is_local) {
1386 CHKERR ISGetSize(is_local, &is_nb_rows);
1387#ifndef NDEBUG
1388 if constexpr (debug_schur) {
1389 CHKERR ISView(is_local, PETSC_VIEWER_STDOUT_WORLD);
1390 }
1391#endif
1392 CHKERR ISGetIndices(is_local, &ptr);
1393 }
1394
1395 int loc_m, loc_n;
1396 CHKERR MatGetLocalSize(A, &loc_m, &loc_n);
1397
1398 for (auto n = 0; n != is_nb_rows; ++n) {
1399 auto row = ptr[n];
1400 auto rlo = view.get<0>().lower_bound(row);
1401 auto rhi = view.get<0>().end();
1402 for (; rlo != rhi; ++rlo) {
1403 auto r_shift = row - (*rlo)->getRow();
1404 if (r_shift >= 0 && r_shift < (*rlo)->getNbRows()) {
1405 auto *ptr = &(*ctx->dataBlocksPtr)[(*rlo)->getMatShift()];
1406 for (auto i = 0; i != (*rlo)->getNbCols(); ++i) {
1407 ptr[i + r_shift * (*rlo)->getNbCols()] = 0;
1408 }
1409 } else if ((*rlo)->getRow() + (*rlo)->getNbRows() > row) {
1410 break;
1411 }
1412 }
1413 }
1414
1415 for (auto n = 0; n != is_nb_rows; ++n) {
1416 auto col = ptr[n];
1417 auto clo = view.get<1>().lower_bound(col);
1418 auto chi = view.get<1>().end();
1419 for (; clo != chi; ++clo) {
1420 auto c_shift = col - (*clo)->getCol();
1421 if (c_shift >= 0 && c_shift < (*clo)->getNbCols()) {
1422
1423 auto *ptr = &(*ctx->dataBlocksPtr)[(*clo)->getMatShift()];
1424 for (auto i = 0; i != (*clo)->getNbRows(); ++i) {
1425 ptr[c_shift + i * (*clo)->getNbCols()] = 0;
1426 }
1427
1428 // diagonal
1429 if (
1430
1431 (*clo)->getRow() == (*clo)->getCol() &&
1432 (*clo)->getLocRow() < loc_m && (*clo)->getLocCol() < loc_n
1433
1434 ) {
1435 auto r_shift = col - (*clo)->getCol();
1436 if (r_shift >= 0 && r_shift < (*clo)->getNbRows()) {
1437 ptr[c_shift + r_shift * (*clo)->getNbCols()] = diag;
1438 }
1439 }
1440 } else if ((*clo)->getCol() + (*clo)->getNbCols() > col) {
1441 break;
1442 }
1443 }
1444 }
1445
1446 if (is_local) {
1447 CHKERR ISRestoreIndices(is_local, &ptr);
1448 }
1449
1450 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_zeroRowsAndCols, 0, 0, 0, 0);
1451
1453}
1454
1455static PetscErrorCode mat_zero(Mat m) {
1457 BlockStructure *ctx;
1458 CHKERR MatShellGetContext(m, (void **)&ctx);
1459 if (ctx->dataBlocksPtr)
1460 std::fill(ctx->dataBlocksPtr->begin(), ctx->dataBlocksPtr->end(), 0.);
1461 if (ctx->preconditionerBlocksPtr)
1462 std::fill(ctx->preconditionerBlocksPtr->begin(),
1463 ctx->preconditionerBlocksPtr->end(), 0.);
1465}
1466
1469 CHKERR MatShellSetManageScalingShifts(mat_raw);
1470 CHKERR MatShellSetOperation(mat_raw, MATOP_MULT, (void (*)(void))mult);
1471 CHKERR MatShellSetOperation(mat_raw, MATOP_MULT_ADD,
1472 (void (*)(void))mult_add);
1473 CHKERR MatShellSetOperation(mat_raw, MATOP_MULT_TRANSPOSE,
1474 (void (*)(void))mult_transpose);
1475 CHKERR MatShellSetOperation(mat_raw, MATOP_MULT_TRANSPOSE_ADD,
1476 (void (*)(void))mult_transpose_add);
1477 CHKERR MatShellSetOperation(mat_raw, MATOP_SOLVE, (void (*)(void))solve);
1478 CHKERR MatShellSetOperation(mat_raw, MATOP_SOLVE_ADD,
1479 (void (*)(void))solve_add);
1480 CHKERR MatShellSetOperation(mat_raw, MATOP_SOLVE_TRANSPOSE,
1481 (void (*)(void))solve_transpose);
1482 CHKERR MatShellSetOperation(mat_raw, MATOP_SOLVE_TRANSPOSE_ADD,
1483 (void (*)(void))solve_transpose_add);
1484 CHKERR MatShellSetOperation(mat_raw, MATOP_ZERO_ENTRIES,
1485 (void (*)(void))mat_zero);
1486 CHKERR MatShellSetOperation(mat_raw, MATOP_ZERO_ROWS_COLUMNS,
1487 (void (*)(void))zero_rows_columns);
1488
1490};
1491
1493 boost::shared_ptr<BlockStructure> data) {
1494
1495 auto problem_ptr = getProblemPtr(dm);
1496 auto nb_local = problem_ptr->nbLocDofsRow;
1497 auto nb_global = problem_ptr->nbDofsRow;
1498
1499 // check in nb, rows is equal to nb. columns.
1500 if (nb_local != problem_ptr->nbLocDofsCol) {
1501 MOFEM_LOG("SELF", Sev::error)
1502 << "Wrong size " << nb_local << " != " << problem_ptr->nbLocDofsCol;
1504 "nb. cols is inconsistent with nb. rows");
1505 }
1506 if (nb_global != problem_ptr->nbDofsCol) {
1507 MOFEM_LOG("SELF", Sev::error)
1508 << "Wrong size " << nb_global << " != " << problem_ptr->nbDofsCol;
1510 "nb. cols is inconsistent with nb. rows");
1511 }
1512
1513 // get comm from DM
1514 MPI_Comm comm;
1515 CHKERR PetscObjectGetComm((PetscObject)dm, &comm);
1516
1517 Mat mat_raw;
1518 CHKERR MatCreateShell(comm, nb_local, nb_local, nb_global, nb_global,
1519 (void *)data.get(), &mat_raw);
1520 CHKERR setSchurBlockMatOps(mat_raw);
1521 // CHKERR PetscObjectSetName((PetscObject)mat_raw, MoFEM_BLOCK_MAT);
1522
1523 return std::make_pair(SmartPetscObj<Mat>(mat_raw), data);
1524}
1525
1527 Mat mat, Vec x, Vec y, InsertMode iora,
1528
1529 boost::function<
1530 int(DiagBlockIndex::BlockIndex::nth_index<0>::type::iterator)>
1531 shift_extractor,
1532
1533 boost::shared_ptr<std::vector<double>> data_blocks_ptr,
1534 bool multiply_by_preconditioner, enum CBLAS_TRANSPOSE trans_A) {
1536 BlockStructure *ctx;
1537 CHKERR MatShellGetContext(mat, (void **)&ctx);
1538
1539 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_BlockStructureMult, 0, 0, 0, 0);
1540
1541 int x_loc_size;
1542 CHKERR VecGetLocalSize(x, &x_loc_size);
1543 int y_loc_size;
1544 CHKERR VecGetLocalSize(y, &y_loc_size);
1545
1546 auto ghost_vecs = std::make_pair(ctx->ghostX, ctx->ghostY);
1547 switch (trans_A) {
1548 case CblasNoTrans:
1549 break;
1550 case CblasTrans:
1551 case CblasConjTrans:
1552 std::swap(ghost_vecs.first, ghost_vecs.second);
1553 break;
1554 };
1555
1556 Vec ghost_x = ghost_vecs.first;
1557 Vec ghost_y = ghost_vecs.second;
1558
1559 CHKERR VecCopy(x, ghost_x);
1560
1561 double *y_array;
1562 Vec loc_ghost_y;
1563 CHKERR VecGhostGetLocalForm(ghost_y, &loc_ghost_y);
1564 int nb_y;
1565 CHKERR VecGetLocalSize(loc_ghost_y, &nb_y);
1566 CHKERR VecGetArray(loc_ghost_y, &y_array);
1567 for (auto i = 0; i != nb_y; ++i)
1568 y_array[i] = 0.;
1569 CHKERR VecRestoreArray(loc_ghost_y, &y_array);
1570 CHKERR VecGhostRestoreLocalForm(ghost_y, &loc_ghost_y);
1571
1572 auto mult = [&](int low_x, int hi_x, int low_y, int hi_y) {
1574
1575 const double *x_array;
1576 Vec loc_ghost_x;
1577 CHKERR VecGhostGetLocalForm(ghost_x, &loc_ghost_x);
1578 CHKERR VecGetArrayRead(loc_ghost_x, &x_array);
1579
1580 double *y_array;
1581 Vec loc_ghost_y;
1582 CHKERR VecGhostGetLocalForm(ghost_y, &loc_ghost_y);
1583 int nb_y;
1584 CHKERR VecGetLocalSize(loc_ghost_y, &nb_y);
1585 CHKERR VecGetArray(loc_ghost_y, &y_array);
1586
1587 double *block_ptr = &*data_blocks_ptr->begin();
1588 auto it = ctx->blockIndex.get<0>().begin();
1589 auto hi = ctx->blockIndex.get<0>().end();
1590
1591 auto get_loc_row_no_trans = [](auto it) { return it->getLocRow(); };
1592 auto get_loc_col_no_trans = [](auto it) { return it->getLocCol(); };
1593 auto get_nb_rows_no_trans = [](auto it) { return it->getNbRows(); };
1594 auto get_nb_cols_no_trans = [](auto it) { return it->getNbCols(); };
1595 auto get_loc_row_trans = [](auto it) { return it->getLocCol(); };
1596 auto get_loc_col_trans = [](auto it) { return it->getLocRow(); };
1597
1598 std::function<int(decltype(it))> get_loc_row;
1599 std::function<int(decltype(it))> get_loc_col;
1600 std::function<int(decltype(it))> get_nb_rows;
1601 std::function<int(decltype(it))> get_nb_cols;
1602
1603 switch (trans_A) {
1604 case CblasNoTrans:
1605 get_loc_row = get_loc_row_no_trans;
1606 get_loc_col = get_loc_col_no_trans;
1607 break;
1608 case CblasTrans:
1609 case CblasConjTrans:
1610 get_loc_row = get_loc_row_trans;
1611 get_loc_col = get_loc_col_trans;
1612 break;
1613 };
1614
1615 while (it != hi) {
1616
1617 auto loc_row = get_loc_row(it);
1618 auto loc_col = get_loc_col(it);
1619
1620 if (loc_row < low_y || loc_row >= hi_y ||
1621 loc_col < low_x || loc_col >= hi_x) {
1622 ++it;
1623 continue;
1624 }
1625
1626 auto nb_rows = get_nb_rows_no_trans(it);
1627 auto nb_cols = get_nb_cols_no_trans(it);
1628 auto x_ptr = &x_array[loc_col];
1629 auto y_ptr = &y_array[loc_row];
1630 auto ptr = &block_ptr[shift_extractor(it)];
1631 cblas_dgemv(CblasRowMajor, trans_A, nb_rows, nb_cols, 1.0, ptr, nb_cols,
1632 x_ptr, 1, 1.0, y_ptr, 1);
1633 ++it;
1634 }
1635
1636 if (multiply_by_preconditioner && ctx->multiplyByPreconditioner) {
1637
1638 if (!ctx->preconditionerBlocksPtr)
1639 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1640 "No parentBlockStructurePtr");
1641
1642 auto preconditioner_ptr = &*ctx->preconditionerBlocksPtr->begin();
1643
1644 auto it = ctx->blockIndex.get<0>().begin();
1645 auto hi = ctx->blockIndex.get<0>().end();
1646
1647 while (it != hi) {
1648 auto loc_row = get_loc_row(it);
1649 auto loc_col = get_loc_col(it);
1650
1651 if (loc_row < low_y || loc_row >= hi_y ||
1652 loc_col < low_x || loc_col >= hi_x) {
1653 ++it;
1654 continue;
1655 }
1656
1657 if (it->getMatShift() != -1) {
1658 auto nb_rows = get_nb_rows_no_trans(it);
1659 auto nb_cols = get_nb_cols_no_trans(it);
1660 auto x_ptr = &x_array[loc_col];
1661 auto y_ptr = &y_array[loc_row];
1662 auto ptr = &preconditioner_ptr[it->getMatShift()];
1663 cblas_dgemv(CblasRowMajor, trans_A, nb_rows, nb_cols, 1.0, ptr,
1664 nb_cols, x_ptr, 1, 1.0, y_ptr, 1);
1665 }
1666
1667 ++it;
1668 }
1669 }
1670
1671 CHKERR VecRestoreArrayRead(loc_ghost_x, &x_array);
1672 CHKERR VecRestoreArray(loc_ghost_y, &y_array);
1673 CHKERR VecGhostRestoreLocalForm(ghost_x, &loc_ghost_x);
1674 CHKERR VecGhostRestoreLocalForm(ghost_y, &loc_ghost_y);
1676 };
1677
1678 constexpr auto max_int = std::numeric_limits<int>::max();
1679 CHKERR VecGhostUpdateBegin(ghost_x, INSERT_VALUES, SCATTER_FORWARD);
1680 CHKERR mult(0, x_loc_size, 0, max_int);
1681 CHKERR VecGhostUpdateEnd(ghost_x, INSERT_VALUES, SCATTER_FORWARD);
1682 CHKERR mult(x_loc_size, max_int, 0, max_int);
1683
1684 CHKERR VecGhostUpdateBegin(ghost_y, ADD_VALUES, SCATTER_REVERSE);
1685 CHKERR VecGhostUpdateEnd(ghost_y, ADD_VALUES, SCATTER_REVERSE);
1686
1687 switch (iora) {
1688 case INSERT_VALUES:
1689 CHKERR VecCopy(ghost_y, y);
1690 break;
1691 case ADD_VALUES:
1692 CHKERR VecAXPY(y, 1., ghost_y);
1693 break;
1694 default:
1695 CHK_THROW_MESSAGE(MOFEM_NOT_IMPLEMENTED, "Wrong InsertMode");
1696 }
1697
1698#ifndef NDEBUG
1699
1700 auto print_norm = [&](auto msg, auto y) {
1702 PetscReal norm;
1703 CHKERR VecNorm(y, NORM_2, &norm);
1704 int nb_loc_y;
1705 CHKERR VecGetLocalSize(y, &nb_loc_y);
1706 int nb_y;
1707 CHKERR VecGetSize(y, &nb_y);
1708 MOFEM_LOG("WORLD", Sev::noisy)
1709 << msg << " " << nb_y << " " << nb_loc_y << " norm " << norm;
1711 };
1712
1713 switch (iora) {
1714 case INSERT_VALUES:
1715 print_norm("mult_schur_block_shell insert x", x);
1716 print_norm("mult_schur_block_shell insert y", y);
1717 break;
1718 case ADD_VALUES:
1719 print_norm("mult_schur_block_shell add x", x);
1720 print_norm("mult_schur_block_shell add y", y);
1721 break;
1722 default:
1723 CHK_THROW_MESSAGE(MOFEM_NOT_IMPLEMENTED, "Wrong InsertMode");
1724 }
1725
1726#endif // NDEBUG
1727
1728 // PetscLogFlops(xxx)
1729 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_BlockStructureMult, 0, 0, 0, 0);
1730
1732}
1733
1734static MoFEMErrorCode solve_schur_block_shell(Mat mat, Vec y, Vec x,
1735 InsertMode iora,
1736 enum CBLAS_TRANSPOSE trans_A) {
1737 using matrix_range = ublas::matrix_range<MatrixDouble>;
1738 using range = ublas::range;
1740 BlockStructure *ctx;
1741 CHKERR MatShellGetContext(mat, (void **)&ctx);
1742
1743 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_BlockStructureSolve, 0, 0, 0, 0);
1744
1745 if (ctx->multiplyByPreconditioner) {
1746 if (!ctx->preconditionerBlocksPtr)
1747 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1748 "No preconditionerBlocksPtr");
1749 }
1750
1751 if (iora == INSERT_VALUES) {
1752 CHKERR VecZeroEntries(x);
1753 }
1754
1755 double *x_array;
1756 CHKERR VecGetArray(x, &x_array);
1757 const double *y_array;
1758 CHKERR VecGetArrayRead(y, &y_array);
1759
1760 std::vector<const DiagBlockIndex::Indexes *> blocks;
1762 VectorDouble block_y;
1763 VectorInt ipiv;
1764
1765 auto get_id = [](auto bi) { return bi->getRowUId(); };
1766 auto get_block_data = [](auto bi) {
1767 return std::make_tuple(bi->getNbRows(), bi->getNbRows(), bi->getLocRow());
1768 };
1769 auto block_sorter = [](auto a, auto b) {
1770 if (a->getRowUId() == b->getRowUId())
1771 return a->getColUId() < b->getColUId();
1772 else
1773 return (a->getRowUId() < b->getRowUId());
1774 };
1775
1776 auto &block_index = ctx->blockIndex.get<1>();
1777 for (auto it = block_index.begin(); it != block_index.end();) {
1778
1779 // get all blocks on finit element
1780 blocks.clear();
1781 auto last_uid = it->getFEUId();
1782 while (it != block_index.end() && it->getFEUId() == last_uid) {
1783 blocks.push_back(&*it);
1784 ++it;
1785 }
1786
1787 // set local indices
1788 std::map<UId, std::tuple<int /*nb rows*/, int /*nb cols*/,
1789 int /*loc index in vector*/>>
1790 block_indexing; // uid block map
1791 for (auto b : blocks) {
1792 // find all unique row uids (columns will be the same)
1793 if (block_indexing.find(get_id(b)) == block_indexing.end()) {
1794 block_indexing[get_id(b)] = get_block_data(b);
1795 }
1796 }
1797 std::sort(blocks.begin(), blocks.end(), block_sorter);
1798
1799 // set indexes to block
1800 int mat_block_size = 0; // size of block matrix
1801 for (auto &b : block_indexing) {
1802 auto &index = std::get<0>(b.second);
1803 auto block_row_size = index;
1804 index = mat_block_size;
1805 mat_block_size += block_row_size;
1806 }
1807
1808 std::vector<std::tuple<int, int, int, int, int>> block_data;
1809 block_data.reserve(blocks.size());
1810 for (auto s : blocks) {
1811 auto ruid = s->getRowUId();
1812 auto cuid = s->getColUId();
1813 auto &rbi = block_indexing.at(ruid);
1814 auto &cbi = block_indexing.at(cuid);
1815 if (auto shift = s->getMatShift(); shift != -1) {
1816 block_data.push_back(std::make_tuple(
1817
1818 shift,
1819
1820 s->getNbRows(), s->getNbCols(),
1821
1822 get<0>(rbi), get<0>(cbi))
1823
1824 );
1825 }
1826 }
1827
1828 // fill block matrix for element
1829 block_mat.resize(mat_block_size, mat_block_size, false);
1830 block_mat.clear();
1831 for (auto &bd : block_data) {
1832 auto &[shift, nb_rows, nb_cols, ridx, cidx] = bd;
1833 auto ptr = &(*(ctx->dataBlocksPtr))[shift];
1834 for (auto i = 0; i != nb_rows; ++i, ptr += nb_cols) {
1835 auto sub_ptr = &block_mat(ridx + i, cidx);
1836 cblas_dcopy(nb_cols, ptr, 1, sub_ptr, 1);
1837 }
1838 }
1839 if (ctx->multiplyByPreconditioner) {
1840 for (auto &bd : block_data) {
1841 auto &[shift, nb_rows, nb_cols, ridx, cidx] = bd;
1842 auto ptr = &(*(ctx->preconditionerBlocksPtr))[shift];
1843 for (auto i = 0; i != nb_rows; ++i, ptr += nb_cols) {
1844 auto sub_ptr = &block_mat(ridx + i, cidx);
1845 cblas_daxpy(nb_cols, 1., ptr, 1, sub_ptr, 1);
1846 }
1847 }
1848 }
1849 if (trans_A == CblasNoTrans)
1850 block_mat = trans(block_mat);
1851
1852 block_y.resize(mat_block_size, false);
1853 block_y.clear();
1854
1855 for (auto &b : block_indexing) {
1856 // note that index for local block_y starts from 0,
1857 // also is positon for rows, or columns depending on
1858 // trans_A
1859 auto [index, size, loc] = b.second;
1860 auto ptr = &y_array[loc];
1861 cblas_dcopy(size, ptr, 1, &block_y[index], 1);
1862 }
1863
1864 // note: this not exploits symmetry, requires more implementation
1865 constexpr int nrhs = 1;
1866 ipiv.resize(mat_block_size, false);
1867 auto info = lapack_dgesv(mat_block_size, nrhs, &*block_mat.data().begin(),
1868 mat_block_size, &*ipiv.data().begin(),
1869 &*block_y.data().begin(), mat_block_size);
1870 if (info != 0) {
1871 MOFEM_LOG("SELF", Sev::error)
1872 << "error lapack solve dgesv info = " << info;
1873 MOFEM_LOG("SELF", Sev::error) << block_mat;
1874 SETERRQ(PETSC_COMM_SELF, MOFEM_OPERATION_UNSUCCESSFUL,
1875 "error lapack solve dgesv info = %d", info);
1876 }
1877
1878 for (auto &b : block_indexing) {
1879 auto [index, size, loc] = b.second;
1880 auto ptr = &x_array[loc];
1881 cblas_daxpy(size, 1.0, &block_y[index], 1, ptr, 1);
1882 }
1883 }
1884
1885 CHKERR VecRestoreArray(x, &x_array);
1886 CHKERR VecRestoreArrayRead(y, &y_array);
1887
1888 // PetscLogFlops(xxx)
1889 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_BlockStructureSolve, 0, 0, 0, 0);
1890
1892}
1893
1896 std::vector<std::string> fields_name,
1897 std::vector<boost::shared_ptr<Range>> field_ents,
1898 SmartPetscObj<AO> ao) {
1899 using matrix_range = ublas::matrix_range<MatrixDouble>;
1900 using range = ublas::range;
1902 BlockStructure *ctx;
1903 CHKERR MatShellGetContext(B, (void **)&ctx);
1904
1905 constexpr bool debug = false;
1906
1907 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_AssembleSchurMat, 0, 0, 0, 0);
1908
1909 if (ao.use_count() == 0) {
1910 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY, "No AO set");
1911 }
1912
1913 if (ctx->multiplyByPreconditioner) {
1914 if (!ctx->preconditionerBlocksPtr)
1915 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1916 "No preconditionerBlocksPtr");
1917 }
1918
1919 std::vector<std::pair<UId, UId>> a00_fields_uids;
1920 a00_fields_uids.reserve(fields_name.size());
1921 auto it_field_ents = field_ents.begin();
1922 if (fields_name.size() != field_ents.size()) {
1923 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
1924 "fields_name.size() != field_ents.size()");
1925 }
1926 for (auto &f : fields_name) {
1927 int bn = m_field.get_field_bit_number(f);
1928 if (*it_field_ents) {
1929 for (auto p = (*it_field_ents)->pair_begin();
1930 p != (*it_field_ents)->pair_end(); ++p) {
1931 a00_fields_uids.emplace_back(
1932 DofEntity::getLoFieldEntityUId(bn, p->first),
1933 DofEntity::getHiFieldEntityUId(bn, p->second));
1934 }
1935 } else {
1936 a00_fields_uids.emplace_back(FieldEntity::getLoBitNumberUId(bn),
1938 }
1939 ++it_field_ents;
1940 }
1941
1942 std::vector<const DiagBlockIndex::Indexes *> fe_blocks;
1943 std::vector<const DiagBlockIndex::Indexes *> a00_blocks;
1944 std::vector<const DiagBlockIndex::Indexes *> a01_blocks;
1945 std::vector<const DiagBlockIndex::Indexes *> a10_blocks;
1946 std::vector<const DiagBlockIndex::Indexes *> a11_blocks;
1947 MatrixDouble block_mat_a00;
1948 MatrixDouble block_mat_a01;
1949 MatrixDouble block_mat_a10;
1950 MatrixDouble block_mat_a11;
1951 VectorInt ipiv;
1952 std::vector<int> glob_idx;
1953 std::vector<int> row_glob_idx, col_glob_idx;
1954
1955 auto &block_index = ctx->blockIndex.get<1>();
1956 for (auto it = block_index.begin(); it != block_index.end();) {
1957
1958 // get blocks on finit element
1959 fe_blocks.clear();
1960 a00_blocks.clear();
1961 a01_blocks.clear();
1962 a10_blocks.clear();
1963 a11_blocks.clear();
1964
1965 auto last_uid = it->getFEUId();
1966 while (it != block_index.end() && it->getFEUId() == last_uid) {
1967 fe_blocks.push_back(&*it);
1968 ++it;
1969 }
1970
1971 for (auto b_ptr : fe_blocks) {
1972 auto check_id = [&](auto uid) {
1973 for (auto &uid_pair : a00_fields_uids) {
1974 if (uid >= uid_pair.first && uid < uid_pair.second) {
1975 return true;
1976 }
1977 }
1978 return false;
1979 };
1980 auto r00 = check_id(b_ptr->getRowUId());
1981 auto c00 = check_id(b_ptr->getColUId());
1982 if (r00 && c00) {
1983 a00_blocks.push_back(b_ptr);
1984 } else if (r00 && !c00) {
1985 a01_blocks.push_back(b_ptr);
1986 } else if (!r00 && c00) {
1987 a10_blocks.push_back(b_ptr);
1988 } else {
1989 a11_blocks.push_back(b_ptr);
1990 }
1991 }
1992
1993 for (auto b : a11_blocks) {
1994
1995 row_glob_idx.resize(b->getNbRows());
1996 col_glob_idx.resize(b->getNbCols());
1997
1998 std::iota(row_glob_idx.begin(), row_glob_idx.end(), b->getRow());
1999 std::iota(col_glob_idx.begin(), col_glob_idx.end(), b->getCol());
2000 CHKERR AOApplicationToPetsc(ao, col_glob_idx.size(), col_glob_idx.data());
2001 CHKERR AOApplicationToPetsc(ao, row_glob_idx.size(), row_glob_idx.data());
2002
2003 auto shift = b->getMatShift();
2004 if (shift != -1) {
2005 auto ptr = &(*(ctx->dataBlocksPtr))[shift];
2006 CHKERR MatSetValues(S, row_glob_idx.size(), row_glob_idx.data(),
2007 col_glob_idx.size(), col_glob_idx.data(), ptr,
2008 ADD_VALUES);
2009 if (ctx->multiplyByPreconditioner) {
2010 auto ptr = &(*(ctx->preconditionerBlocksPtr))[shift];
2011 CHKERR MatSetValues(S, row_glob_idx.size(), row_glob_idx.data(),
2012 col_glob_idx.size(), col_glob_idx.data(), ptr,
2013 ADD_VALUES);
2014 }
2015 }
2016 }
2017
2018 if (debug) {
2019 MOFEM_LOG("WORLD", Sev::warning)
2020 << "a00_blocks.size() " << a00_blocks.size();
2021 MOFEM_LOG("WORLD", Sev::warning)
2022 << "a01_blocks.size() " << a01_blocks.size();
2023 MOFEM_LOG("WORLD", Sev::warning)
2024 << "a10_blocks.size() " << a10_blocks.size();
2025 MOFEM_LOG("WORLD", Sev::warning)
2026 << "a11_blocks.size() " << a11_blocks.size();
2027 }
2028
2029 if(a00_blocks.size()) {
2030
2031 for (auto r : a00_blocks) {
2032 auto r_uid = r->getRowUId();
2033 auto range = ctx->blockIndex.get<2>().equal_range(r_uid);
2034 for (auto it = range.first; it != range.second; ++it) {
2035 if (it->getFEUId() != last_uid && it->getColUId() != r_uid) {
2036 a01_blocks.push_back(&*it);
2037 }
2038 }
2039 }
2040
2041 for (auto c : a00_blocks) {
2042 auto c_uid = c->getColUId();
2043 auto range = ctx->blockIndex.get<3>().equal_range(c_uid);
2044 for (auto it = range.first; it != range.second; ++it) {
2045 if (it->getFEUId() != last_uid && it->getRowUId() != c_uid) {
2046 a10_blocks.push_back(&*it);
2047 }
2048 }
2049 }
2050
2051 auto sort = [](auto &blocks) {
2052 std::sort(blocks.begin(), blocks.end(), [](auto a, auto b) {
2053 if (a->getRowUId() == b->getRowUId())
2054 return a->getColUId() < b->getColUId();
2055 else
2056 return a->getRowUId() < b->getRowUId();
2057 });
2058 };
2059
2060 sort(a00_blocks);
2061 sort(a01_blocks);
2062 sort(a10_blocks);
2063
2064 if (debug) {
2065 MOFEM_LOG("WORLD", Sev::warning)
2066 << "a01_blocks.size() " << a01_blocks.size();
2067 MOFEM_LOG("WORLD", Sev::warning)
2068 << "a10_blocks.size() " << a10_blocks.size();
2069 }
2070
2071 // set local indices
2072 auto set_local_indices = [](auto &row_blocks, auto &col_blocks) {
2073 // index, size, global index
2074 std::map<UId, std::tuple<int, int, int>>
2075 block_indexing; // uid block map
2076 for (auto b : row_blocks) {
2077 if (block_indexing.find(b->getRowUId()) == block_indexing.end()) {
2078 block_indexing[b->getRowUId()] =
2079 std::make_tuple(b->getNbRows(), b->getNbRows(), b->getRow());
2080 }
2081 }
2082 for (auto b : col_blocks) {
2083 if (block_indexing.find(b->getColUId()) == block_indexing.end()) {
2084 block_indexing[b->getColUId()] =
2085 std::make_tuple(b->getNbCols(), b->getNbCols(), b->getCol());
2086 }
2087 }
2088 // set indexes to block
2089 int mat_block_size = 0; // size of block matrix
2090 for (auto &b : block_indexing) {
2091 // index, size, global index
2092 auto &[index, size, glob] = b.second;
2093 index = mat_block_size;
2094 mat_block_size += size;
2095 }
2096 return std::make_pair(block_indexing, mat_block_size);
2097 };
2098
2099 auto set_block = [&](auto &blocks, auto &row_indexing, auto &col_indexing,
2100 auto &block_mat, int row_block_size,
2101 int col_block_size) {
2102 std::vector<std::tuple<int, int, int, int, int>> block_data;
2103 block_data.reserve(blocks.size());
2104 for (auto s : blocks) {
2105 auto ruid = s->getRowUId();
2106 auto cuid = s->getColUId();
2107 auto &rbi = row_indexing.at(ruid);
2108 auto &cbi = col_indexing.at(cuid);
2109 if (auto shift = s->getMatShift(); shift != -1) {
2110 block_data.push_back(std::make_tuple(
2111
2112 shift,
2113
2114 s->getNbRows(), s->getNbCols(),
2115
2116 get<0>(rbi), get<0>(cbi))
2117
2118 );
2119 }
2120 }
2121 block_mat.resize(row_block_size, col_block_size, false);
2122 block_mat.clear();
2123 for (auto &bd : block_data) {
2124 auto &[shift, nb_rows, nb_cols, ridx, cidx] = bd;
2125 auto ptr = &(*(ctx->dataBlocksPtr))[shift];
2126 for (auto i = 0; i != nb_rows; ++i, ptr += nb_cols) {
2127 auto sub_ptr = &block_mat(ridx + i, cidx);
2128 cblas_dcopy(nb_cols, ptr, 1, sub_ptr, 1);
2129 }
2130 }
2131 if (ctx->multiplyByPreconditioner) {
2132 for (auto &bd : block_data) {
2133 auto &[shift, nb_rows, nb_cols, ridx, cidx] = bd;
2134 auto ptr = &(*(ctx->preconditionerBlocksPtr))[shift];
2135 for (auto i = 0; i != nb_rows; ++i, ptr += nb_cols) {
2136 auto sub_ptr = &block_mat(ridx + i, cidx);
2137 cblas_daxpy(nb_cols, 1., ptr, 1, sub_ptr, 1);
2138 }
2139 }
2140 }
2141 };
2142
2143 auto [a00_indexing, a00_size] = set_local_indices(a00_blocks, a00_blocks);
2144 auto [a11_indexing, a11_size] = set_local_indices(a10_blocks, a01_blocks);
2145
2146 if (debug) {
2147 MOFEM_LOG("WORLD", Sev::warning)
2148 << "a00_indexing.size() " << a00_indexing.size() << " a00_size "
2149 << a00_size;
2150 MOFEM_LOG("WORLD", Sev::warning)
2151 << "a11_indexing.size() " << a11_indexing.size() << " a11_size "
2152 << a11_size;
2153 }
2154
2155 set_block(a00_blocks, a00_indexing, a00_indexing, block_mat_a00, a00_size,
2156 a00_size);
2157 set_block(a01_blocks, a00_indexing, a11_indexing, block_mat_a01, a00_size,
2158 a11_size);
2159 set_block(a10_blocks, a11_indexing, a00_indexing, block_mat_a10, a11_size,
2160 a00_size);
2161 block_mat_a11.resize(a11_size, a11_size, false);
2162 block_mat_a11.clear();
2163
2164 block_mat_a00 = trans(block_mat_a00);
2165 block_mat_a01 = trans(block_mat_a01);
2166
2167 ipiv.resize(block_mat_a00.size1(), false);
2168 auto info =
2169 lapack_dgesv(block_mat_a00.size1(), block_mat_a01.size1(),
2170 block_mat_a00.data().data(), block_mat_a00.size2(),
2171 &*ipiv.data().begin(), block_mat_a01.data().data(),
2172 block_mat_a01.size2());
2173 if (info)
2174 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2175 "lapack error info = %d", info);
2176
2177 // block_mat_a01 = trans(block_mat_a01);
2178 cblas_dgemm(CblasRowMajor, CblasNoTrans, CblasTrans,
2179 block_mat_a11.size1(), block_mat_a11.size2(),
2180 block_mat_a00.size1(), -1.,
2181
2182 block_mat_a10.data().data(), block_mat_a10.size2(),
2183
2184 block_mat_a01.data().data(), block_mat_a01.size2(),
2185
2186 0., block_mat_a11.data().data(), block_mat_a11.size2()
2187
2188 );
2189
2190 int idx = 0;
2191 std::vector<int> glob_idx(block_mat_a11.size1());
2192 for (auto &r : a11_indexing) {
2193 auto [r_index, r_size, r_glob] = r.second;
2194 std::iota(&glob_idx[idx], &glob_idx[idx + r_size], r_glob);
2195 idx += r_size;
2196 }
2197 CHKERR AOApplicationToPetsc(ao, glob_idx.size(), glob_idx.data());
2198 CHKERR MatSetValues(S, glob_idx.size(), glob_idx.data(), glob_idx.size(),
2199 glob_idx.data(), block_mat_a11.data().data(),
2200 ADD_VALUES);
2201 }
2202 }
2203
2204 // PetscLogFlops(xxx)
2205 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_AssembleSchurMat, 0, 0, 0, 0);
2206
2208}
2209
2210boost::shared_ptr<std::vector<double>> getBlockMatStorageMat(Mat B) {
2211 BlockStructure *ctx;
2212 CHKERR MatShellGetContext(B, (void **)&ctx);
2213 return ctx->dataBlocksPtr;
2214}
2215
2216boost::shared_ptr<std::vector<double>>
2218 BlockStructure *ctx;
2219 CHKERR MatShellGetContext(B, (void **)&ctx);
2220 return ctx->dataBlocksPtr;
2221}
2222
2224 BlockStructure *ctx, const EntitiesFieldData::EntData &row_data,
2225 const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat,
2226 InsertMode iora,
2227 boost::function<int(const DiagBlockIndex::Indexes *)> shift_extractor,
2228 boost::shared_ptr<std::vector<double>> data_blocks_ptr) {
2229
2231
2232 if (row_data.getIndices().empty())
2234 if (col_data.getIndices().empty())
2236
2237#ifndef NDEBUG
2238
2239 PetscLogEventBegin(SchurEvents::MOFEM_EVENT_BlockStructureSetValues, 0, 0, 0,
2240 0);
2241
2242#endif // NDEBUG
2243
2244 switch (iora) {
2245 case ADD_VALUES:
2246 case INSERT_VALUES:
2247 break;
2248 default:
2249 SETERRQ(PETSC_COMM_SELF, MOFEM_NOT_IMPLEMENTED, "Wrong InsertMode");
2250 }
2251
2252 auto get_row_data = [&]() -> std::pair<bool, VectorInt> {
2253 if (auto e_ptr = row_data.getFieldEntities()[0]) {
2254 if (auto stored_data_ptr =
2255 e_ptr->getSharedStoragePtr<EssentialBcStorage>()) {
2256 MOFEM_LOG("SELF", Sev::warning) << "Can lead to unhandled behaviour";
2257 return std::make_pair(true, stored_data_ptr->entityIndices);
2258 }
2259 }
2260 return std::make_pair(false, row_data.getIndices());
2261 };
2262
2263 auto row_indices = get_row_data();
2264
2265 std::vector<int> ent_row_indices;
2266 std::vector<int> ent_col_indices;
2267
2268 for (auto &rent : row_data.getFieldEntities()) {
2269 if (auto r_cache = rent->entityCacheRowDofs.lock()) {
2270
2271 auto &row_uid = rent->getLocalUniqueId();
2272 auto &row_ind = row_indices.second;
2273
2274 for (auto &cent : col_data.getFieldEntities()) {
2275 if (auto c_cache = cent->entityCacheColDofs.lock()) {
2276
2277 auto &col_uid = cent->getLocalUniqueId();
2278 auto &col_ind = col_data.getIndices();
2279
2280 auto it = ctx->blockIndex.get<0>().find(
2281
2282 boost::make_tuple(row_uid, col_uid)
2283
2284 );
2285
2286#ifndef NDEBUG
2287
2288 if (it == ctx->blockIndex.get<0>().end()) {
2289 MOFEM_LOG_CHANNEL("SELF");
2290 MOFEM_TAG_AND_LOG("SELF", Sev::error, "BlockMat")
2291 << "missing block: row "
2292 << row_data.getFieldDofs()[0]->getName() << " col "
2293 << col_data.getFieldDofs()[0]->getName();
2294 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2295 "Block not allocated");
2296 }
2297
2298#endif
2299
2300 auto ri = row_ind.begin();
2301 auto rie = row_ind.end();
2302
2303 auto shift = shift_extractor(&*it);
2304 auto s_mat = &(*data_blocks_ptr)[shift];
2305
2306 auto get_ent_indices = [](auto &cache, auto &ind) {
2307 ind.clear();
2308 ind.reserve(std::distance(cache->loHi[0], cache->loHi[1]));
2309 for (auto e = cache->loHi[0]; e != cache->loHi[1]; ++e) {
2310 auto glob = (*e)->getPetscGlobalDofIdx();
2311 if (glob != -1)
2312 ind.push_back(glob);
2313 }
2314 };
2315
2316 get_ent_indices(r_cache, ent_row_indices);
2317 if (ent_row_indices.empty())
2318 continue;
2319 get_ent_indices(c_cache, ent_col_indices);
2320 if (ent_col_indices.empty())
2321 continue;
2322
2323 if (mat.size1() == ent_row_indices.size() &&
2324 mat.size2() == ent_col_indices.size()) {
2325
2326 if (iora == ADD_VALUES) {
2327 cblas_daxpy(mat.data().size(), 1.0, mat.data().data(), 1, s_mat,
2328 1);
2329 } else {
2330 cblas_dcopy(mat.data().size(), mat.data().data(), 1, s_mat, 1);
2331 }
2332
2333 } else {
2334
2335 int row = 0;
2336 for (auto re : ent_row_indices) {
2337 ri = std::find(ri, rie, re);
2338 if (!(ri == rie && *ri != -1)) {
2339
2340 auto ci = col_ind.begin();
2341 auto cie = col_ind.end();
2342
2343 int col = 0;
2344 for (auto ce : ent_col_indices) {
2345 ci = std::find(ci, cie, ce);
2346 if (!(ci == cie && *ci != -1)) {
2347 auto &m = s_mat[row * ent_col_indices.size() + col];
2348 if (iora == ADD_VALUES) {
2349 m += mat(std::distance(row_ind.begin(), ri),
2350 std::distance(col_ind.begin(), ci));
2351 } else {
2352 m = mat(std::distance(row_ind.begin(), ri),
2353 std::distance(col_ind.begin(), ci));
2354 }
2355 }
2356 ++col;
2357 } // cols
2358 }
2359 ++row;
2360 } // rows
2361 }
2362
2363 } // iterate entities
2364 }
2365 }
2366 }
2367
2368#ifndef NDEBUG
2369 PetscLogEventEnd(SchurEvents::MOFEM_EVENT_BlockStructureSetValues, 0, 0, 0,
2370 0);
2371#endif // NDEBUG
2372
2374}
2375
2378 const EntitiesFieldData::EntData &col_data,
2379 const MatrixDouble &mat, InsertMode iora) {
2381 PetscBool is_mat_shell = PETSC_FALSE;
2382 PetscObjectTypeCompare((PetscObject)M, MATSHELL, &is_mat_shell);
2383 if (is_mat_shell) {
2384 BlockStructure *ctx;
2385 CHKERR MatShellGetContext(M, (void **)&ctx);
2387 ctx, row_data, col_data, mat, iora,
2388 [](const DiagBlockIndex::Indexes *idx) { return idx->getMatShift(); },
2389 ctx->dataBlocksPtr);
2390 } else {
2391 CHKERR MatSetValues<AssemblyTypeSelector<PETSC>>(M, row_data, col_data, mat,
2392 iora);
2393 }
2395}
2396
2398 Mat M, const EntitiesFieldData::EntData &row_data,
2399 const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat,
2400 InsertMode iora) {
2402 PetscBool is_mat_shell = PETSC_FALSE;
2403 PetscObjectTypeCompare((PetscObject)M, MATSHELL, &is_mat_shell);
2404 if (is_mat_shell) {
2405 BlockStructure *ctx;
2406 CHKERR MatShellGetContext(M, (void **)&ctx);
2407 if (!ctx->preconditionerBlocksPtr)
2408 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2409 "No preconditionerBlocksPtr");
2411 ctx, row_data, col_data, mat, iora,
2412 [](const DiagBlockIndex::Indexes *idx) { return idx->getMatShift(); },
2414 } else {
2415 CHKERR MatSetValues<AssemblyTypeSelector<PETSC>>(M, row_data, col_data, mat,
2416 iora);
2417 }
2419}
2420
2421boost::shared_ptr<NestSchurData> createSchurNestedMatrixStruture(
2422
2423 std::pair<SmartPetscObj<DM>, SmartPetscObj<DM>> dms,
2424 boost::shared_ptr<BlockStructure> block_mat_data_ptr,
2425
2426 std::vector<std::string> fields_names, //< a00 fields
2427 std::vector<boost::shared_ptr<Range>>
2428 field_ents, //< a00 ranges (can be null),
2429 bool add_preconditioner_block) {
2430
2431 if (!block_mat_data_ptr) {
2432 CHK_THROW_MESSAGE(MOFEM_DATA_INCONSISTENCY, "Block data not set");
2433 }
2434
2435 if (fields_names.size() != field_ents.size())
2437 "fields_names.size() != field_ents.size()");
2438
2439 auto [schur_dm, block_dm] = dms;
2440 auto schur_prb = getProblemPtr(schur_dm);
2441 auto block_prb = getProblemPtr(block_dm);
2442 // auto m_field_ptr = getInterfacePtr(block_dm);
2443
2444 auto schur_dofs_row = schur_prb->getNumeredRowDofsPtr();
2445 auto schur_dofs_col = schur_prb->getNumeredColDofsPtr();
2446 auto block_dofs_row = block_prb->getNumeredRowDofsPtr();
2447 auto block_dofs_col = block_prb->getNumeredColDofsPtr();
2448
2449 auto ao_schur_row = schur_prb->getSubData()->getSmartRowMap();
2450 auto ao_schur_col = schur_prb->getSubData()->getSmartColMap();
2451 auto ao_block_row = block_prb->getSubData()->getSmartRowMap();
2452 auto ao_block_col = block_prb->getSubData()->getSmartColMap();
2453
2454 auto schur_vec_x = createDMVector(schur_dm);
2455 auto block_vec_x = createDMVector(block_dm);
2456 auto schur_vec_y = vectorDuplicate(schur_vec_x);
2457 auto block_vec_y = vectorDuplicate(block_vec_x);
2458 CHKERR VecSetDM(schur_vec_x, PETSC_NULLPTR);
2459 CHKERR VecSetDM(block_vec_x, PETSC_NULLPTR);
2460 CHKERR VecSetDM(schur_vec_y, PETSC_NULLPTR);
2461 CHKERR VecSetDM(block_vec_y, PETSC_NULLPTR);
2462
2463 auto find_field_ent = [&](auto uid, auto prb, auto rc) {
2464 boost::shared_ptr<NumeredDofEntity_multiIndex> dofs;
2465
2466 switch (rc) {
2467 case ROW:
2468 dofs = prb->getNumeredRowDofsPtr();
2469 break;
2470 case COL:
2471 dofs = prb->getNumeredColDofsPtr();
2472 break;
2473 default:
2475 break;
2476 }
2477
2478 auto lo = dofs->get<Unique_mi_tag>().lower_bound(uid);
2479 if (lo == dofs->get<Unique_mi_tag>().end())
2480 return boost::shared_ptr<NumeredDofEntity>();
2481 auto hi = dofs->get<Unique_mi_tag>().upper_bound(
2483 if (lo != hi)
2484 return *lo;
2485
2486 return boost::shared_ptr<NumeredDofEntity>();
2487 };
2488
2489 std::array<boost::shared_ptr<BlockStructure>, 4> data_ptrs;
2490
2491 for (auto r = 0; r != 3; ++r) {
2492 data_ptrs[r] = boost::make_shared<BlockStructure>();
2493 data_ptrs[r]->dataBlocksPtr = block_mat_data_ptr->dataBlocksPtr;
2494 }
2495 data_ptrs[3] = boost::make_shared<BlockStructure>();
2496 data_ptrs[3]->dataBlocksPtr = block_mat_data_ptr->dataBlocksPtr;
2497
2498 data_ptrs[0]->ghostX = schur_vec_x;
2499 data_ptrs[0]->ghostY = schur_vec_y;
2500 data_ptrs[1]->ghostX = block_vec_x;
2501 data_ptrs[1]->ghostY = schur_vec_y;
2502 data_ptrs[2]->ghostX = schur_vec_x;
2503 data_ptrs[2]->ghostY = block_vec_y;
2504 data_ptrs[3]->ghostX = block_vec_x;
2505 data_ptrs[3]->ghostY = block_vec_y;
2506
2507 int idx = 0;
2508 for (auto &d : block_mat_data_ptr->blockIndex.get<0>()) {
2509
2510 auto insert = [&](auto &m, auto &dof_r, auto &dof_c, auto &d) {
2511 m.insert(
2512
2514 d.getRowUId(), d.getColUId(), d.getFEUId(),
2515
2516 dof_r->getPetscGlobalDofIdx(), dof_c->getPetscGlobalDofIdx(),
2517
2518 d.getNbRows(), d.getNbCols(),
2519
2520 dof_r->getPetscLocalDofIdx(), dof_c->getPetscLocalDofIdx(),
2521
2522 d.getMatShift()}
2523
2524 );
2525 };
2526
2527 auto dof_schur_row = find_field_ent(d.getRowUId(), schur_prb, ROW);
2528 auto dof_schur_col = find_field_ent(d.getColUId(), schur_prb, COL);
2529 auto dof_block_row = find_field_ent(d.getRowUId(), block_prb, ROW);
2530 auto dof_block_col = find_field_ent(d.getColUId(), block_prb, COL);
2531
2532 if (dof_schur_row && dof_schur_col) {
2533 insert(data_ptrs[0]->blockIndex, dof_schur_row, dof_schur_col, d);
2534 }
2535
2536 if (dof_schur_row && dof_block_col) {
2537 insert(data_ptrs[1]->blockIndex, dof_schur_row, dof_block_col, d);
2538 }
2539
2540 if (dof_block_row && dof_schur_col) {
2541 insert(data_ptrs[2]->blockIndex, dof_block_row, dof_schur_col, d);
2542 }
2543
2544 if (dof_block_row && dof_block_col) {
2545 insert(data_ptrs[3]->blockIndex, dof_block_row, dof_block_col, d);
2546 }
2547
2548 ++idx;
2549 }
2550
2551 // set data for a00 solve (inverse blocks)
2552 auto set_up_a00_data = [&](auto inv_block_data) {
2554
2555 if (add_preconditioner_block) {
2556 auto preconditioned_block = boost::make_shared<std::vector<double>>(
2557 inv_block_data->dataBlocksPtr->size(), 0);
2558 inv_block_data->preconditionerBlocksPtr = preconditioned_block;
2559 inv_block_data->multiplyByPreconditioner = true;
2560 block_mat_data_ptr->preconditionerBlocksPtr =
2561 inv_block_data->preconditionerBlocksPtr;
2562 block_mat_data_ptr->multiplyByPreconditioner = false;
2563 }
2564
2566 };
2567
2568 CHKERR set_up_a00_data(data_ptrs[3]);
2569
2570 MPI_Comm comm;
2571 CHKERR PetscObjectGetComm((PetscObject)schur_dm, &comm);
2572
2573 auto create_shell_mat = [&](auto nb_r_loc, auto nb_c_loc, auto nb_r_glob,
2574 auto nb_c_glob, auto data_ptr) {
2575 Mat mat_raw;
2576 CHKERR MatCreateShell(comm, nb_r_loc, nb_c_loc, nb_r_glob, nb_c_glob,
2577 (void *)data_ptr.get(), &mat_raw);
2578 CHKERR setSchurBlockMatOps(mat_raw);
2579 return SmartPetscObj<Mat>(mat_raw);
2580 };
2581
2582 auto schur_nb_global = schur_prb->getNbDofsRow();
2583 auto block_nb_global = block_prb->getNbDofsRow();
2584 auto schur_nb_local = schur_prb->getNbLocalDofsRow();
2585 auto block_nb_local = block_prb->getNbLocalDofsRow();
2586
2587 std::array<SmartPetscObj<Mat>, 4> mats_array;
2588 mats_array[0] =
2589 create_shell_mat(schur_nb_local, schur_nb_local, schur_nb_global,
2590 schur_nb_global, data_ptrs[0]);
2591 mats_array[1] =
2592 create_shell_mat(schur_nb_local, block_nb_local, schur_nb_global,
2593 block_nb_global, data_ptrs[1]);
2594 mats_array[2] =
2595 create_shell_mat(block_nb_local, schur_nb_local, block_nb_global,
2596 schur_nb_global, data_ptrs[2]);
2597 mats_array[3] =
2598 create_shell_mat(block_nb_local, block_nb_local, block_nb_global,
2599 block_nb_global, data_ptrs[3]);
2600
2601 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "NestedSchur")
2602 << "(0, 0) " << schur_nb_local << " " << schur_nb_global << " "
2603 << data_ptrs[0]->blockIndex.size();
2604 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "NestedSchur")
2605 << "(0, 1) " << schur_nb_local << " " << block_nb_local << " "
2606 << schur_nb_global << " " << block_nb_global << " "
2607 << data_ptrs[1]->blockIndex.size();
2608 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "NestedSchur")
2609 << "(1, 0) " << block_nb_local << " " << schur_nb_local << " "
2610 << block_nb_global << " " << schur_nb_global << " "
2611 << data_ptrs[2]->blockIndex.size();
2612 MOFEM_TAG_AND_LOG("SYNC", Sev::verbose, "NestedSchur")
2613 << "(1, 1) " << block_nb_local << " " << block_nb_global << " "
2614 << data_ptrs[3]->blockIndex.size();
2615
2616 MOFEM_LOG_SEVERITY_SYNC(comm, Sev::verbose);
2617
2618 auto schur_is = schur_prb->getSubData()->getSmartRowIs();
2619 auto block_is = block_prb->getSubData()->getSmartRowIs();
2620
2621 return boost::make_shared<NestSchurData>(
2622
2623 mats_array, data_ptrs, block_mat_data_ptr,
2624 std::make_pair(schur_is, block_is)
2625
2626 );
2627}
2628
2629std::pair<SmartPetscObj<Mat>, boost::shared_ptr<NestSchurData>>
2630createSchurNestedMatrix(boost::shared_ptr<NestSchurData> schur_net_data_ptr) {
2631
2632 if (!schur_net_data_ptr)
2634
2635 auto [mat_arrays, data_ptrs, block_mat_data_ptr, iss] = *schur_net_data_ptr;
2636 auto [schur_is, block_is] = iss;
2637
2638 std::array<IS, 2> is_a = {schur_is, block_is};
2639 std::array<Mat, 4> mats_a = {
2640
2641 mat_arrays[0], mat_arrays[1], mat_arrays[2], mat_arrays[3]
2642
2643 };
2644
2645 MPI_Comm comm;
2646 CHKERR PetscObjectGetComm((PetscObject)mat_arrays[0], &comm);
2647
2648 Mat mat_raw;
2649 CHKERR MatCreateNest(
2650
2651 comm, 2, is_a.data(), 2, is_a.data(), mats_a.data(), &mat_raw
2652
2653 );
2654
2655 return std::make_pair(SmartPetscObj<Mat>(mat_raw), schur_net_data_ptr);
2656}
2657
2661
2662OpSchurAssembleBase *
2663createOpSchurAssembleEnd(std::vector<std::string> fields_name,
2664 std::vector<boost::shared_ptr<Range>> field_ents,
2666 bool sym_schur, bool symm_op) {
2667 if (symm_op)
2668 return new OpSchurAssembleEnd<SchurDSYSV>(fields_name, field_ents, ao,
2669 schur, sym_schur, symm_op);
2670 else
2671 return new OpSchurAssembleEnd<SchurDGESV>(fields_name, field_ents, ao,
2672 schur, sym_schur, symm_op);
2673}
2674
2676 boost::shared_ptr<std::vector<unsigned char>> marker_ptr, double diag_val) {
2677 return new OpSchurZeroRowsAndCols(marker_ptr, diag_val);
2678}
2679
2680OpSchurAssembleBase *
2681createOpSchurAssembleEnd(std::vector<std::string> fields_name,
2682 std::vector<boost::shared_ptr<Range>> field_ents,
2683 std::vector<SmartPetscObj<AO>> sequence_of_aos,
2684 std::vector<SmartPetscObj<Mat>> sequence_of_mats,
2685 std::vector<bool> sym_schur, bool symm_op,
2686 boost::shared_ptr<BlockStructure> diag_blocks) {
2688 fields_name, field_ents, sequence_of_aos.back(), sequence_of_mats.back(),
2689 sym_schur.back(), symm_op);
2690}
2691
2692OpSchurAssembleBase *
2693createOpSchurAssembleEnd(std::vector<std::string> fields_name,
2694 std::vector<boost::shared_ptr<Range>> field_ents,
2695 std::vector<SmartPetscObj<AO>> sequence_of_aos,
2696 std::vector<SmartPetscObj<Mat>> sequence_of_mats,
2697 std::vector<bool> sym_schur,
2698 std::vector<double> diag_eps, bool symm_op,
2699 boost::shared_ptr<BlockStructure> diag_blocks) {
2701 fields_name, field_ents, sequence_of_aos.back(), sequence_of_mats.back(),
2702 sym_schur.back(), symm_op);
2703}
2704
2707
2708 auto apply = [](PC pc, Vec f, Vec x) {
2710 Mat A, P;
2711 CHKERR PCGetOperators(pc, &A, &P);
2712 CHKERR MatSolve(P, f, x);
2714 };
2715
2716 auto apply_transpose = [](PC pc, Vec f, Vec x) {
2718 Mat A, P;
2719 CHKERR PCGetOperators(pc, &A, &P);
2720 CHKERR MatSolveTranspose(P, f, x);
2722 };
2723
2724 CHKERR PCSetType(pc, PCSHELL);
2725 CHKERR PCShellSetApply(pc, apply);
2726 CHKERR PCShellSetApplyTranspose(pc, apply_transpose);
2727 CHKERR PCShellSetName(pc, "MoFEMSchurBlockPC");
2728
2730}
2731
2732// This is used to assemble local element matrices for Schur complement
2733// and matrix for KSP
2734template <>
2737 const EntitiesFieldData::EntData &col_data,
2738 const MatrixDouble &mat, InsertMode iora) {
2739 return SchurElemMats::MatSetValues(M, row_data, col_data, mat, iora);
2740}
2741
2744 const EntitiesFieldData::EntData &col_data,
2745 const MatrixDouble &mat, InsertMode iora) {
2746 return MatSetValues<AssemblyTypeSelector<PETSC>>(M, row_data, col_data, mat,
2747 iora);
2748}
2749
2750// Is assumed that standard PETSc assembly works for matrices used by KSP
2753
2754// We assemble matrix for KSP and store local matrices for Schur complement.
2755// Schur complement is calculated and assembled in OpSchurAssembleEnd.
2758 const EntitiesFieldData::EntData &col_data,
2759 const MatrixDouble &mat, InsertMode iora) {
2761 CHKERR assembleStorage(row_data, col_data, mat, iora);
2763 mat, iora);
2765}
2766
2767// All is now wrapped in specialisation of
2768// MatSetValues<AssemblyTypeSelector<SCHUR>>
2769template <>
2770MoFEMErrorCode MatSetValues<AssemblyTypeSelector<SCHUR>>(
2771 Mat M, const EntitiesFieldData::EntData &row_data,
2772 const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat,
2773 InsertMode iora) {
2774 return MatSetValues<SchurElemMats>(M, row_data, col_data, mat, iora);
2775}
2776
2777// Standard (PETSc) assembly block matrices do noe work
2778template <>
2781 const EntitiesFieldData::EntData &col_data,
2782 const MatrixDouble &mat, InsertMode iora) {
2783 return shell_block_mat_asmb_wrap(M, row_data, col_data, mat, iora);
2784}
2785
2787
2788 static MoFEMErrorCode MatSetValues(Mat M,
2789 const EntitiesFieldData::EntData &row_data,
2790 const EntitiesFieldData::EntData &col_data,
2791 const MatrixDouble &mat, InsertMode iora);
2792};
2793
2797
2800 const EntitiesFieldData::EntData &row_data,
2801 const EntitiesFieldData::EntData &col_data,
2802 const MatrixDouble &mat, InsertMode iora) {
2804 CHKERR assembleStorage(row_data, col_data, mat, iora);
2806 col_data, mat, iora);
2808}
2809
2810template <>
2813 const EntitiesFieldData::EntData &row_data,
2814 const EntitiesFieldData::EntData &col_data,
2815 const MatrixDouble &mat, InsertMode iora) {
2816 return SchurElemMatsBlock::MatSetValues(M, row_data, col_data, mat, iora);
2817}
2818
2820
2821 static MoFEMErrorCode MatSetValues(Mat M,
2822 const EntitiesFieldData::EntData &row_data,
2823 const EntitiesFieldData::EntData &col_data,
2824 const MatrixDouble &mat, InsertMode iora);
2825};
2826
2830
2832 Mat M, const EntitiesFieldData::EntData &row_data,
2833 const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat,
2834 InsertMode iora) {
2836 CHKERR assembleStorage(row_data, col_data, mat, iora);
2838 M, row_data, col_data, mat, iora);
2840}
2841
2842template <>
2844 Mat M, const EntitiesFieldData::EntData &row_data,
2845 const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat,
2846 InsertMode iora) {
2847 return SchurElemMatsPreconditionedBlock::MatSetValues(M, row_data, col_data,
2848 mat, iora);
2849}
2850
2852schurSwitchPreconditioner(boost::shared_ptr<BlockStructure> block_mat_data) {
2854 if (block_mat_data->multiplyByPreconditioner) {
2855 block_mat_data->multiplyByPreconditioner = false;
2856 } else {
2857 if (!block_mat_data->preconditionerBlocksPtr) {
2858 SETERRQ(PETSC_COMM_SELF, MOFEM_DATA_INCONSISTENCY,
2859 "preconditionerBlocksPtr not set");
2860 }
2861 block_mat_data->multiplyByPreconditioner = true;
2862 }
2864}
2865
2867schurSaveBlockMesh(boost::shared_ptr<BlockStructure> block_mat_data,
2868 std::string filename) {
2870
2871 moab::Core core;
2872 moab::Interface &moab = core;
2873
2874 ReadUtilIface *iface;
2875 CHKERR moab.query_interface(iface);
2876
2877 auto size = block_mat_data->blockIndex.size();
2878
2879 EntityHandle starting_vertex;
2880 vector<double *> arrays_coord;
2881 CHKERR iface->get_node_coords(3, 4 * size, 0, starting_vertex, arrays_coord);
2882 EntityHandle starting_handle;
2883 EntityHandle *array;
2884 CHKERR iface->get_element_connect(size, 4, MBQUAD, 0, starting_handle, array);
2885
2886 double def_val_nrm2 = 0;
2887 Tag th_nrm2;
2888 CHKERR moab.tag_get_handle("nrm2", 1, MB_TYPE_DOUBLE, th_nrm2,
2889 MB_TAG_CREAT | MB_TAG_DENSE, &def_val_nrm2);
2890
2891 int def_val_mat_shift = 0;
2892 Tag th_mat_shift;
2893 CHKERR moab.tag_get_handle("mat_shift", 1, MB_TYPE_INTEGER, th_mat_shift,
2894 MB_TAG_CREAT | MB_TAG_DENSE, &def_val_mat_shift);
2895
2896 int i = 0;
2897 for (auto &d : block_mat_data->blockIndex) {
2898 auto row = d.getRow();
2899 auto col = d.getCol();
2900 auto nb_rows = d.getNbRows();
2901 auto nb_cols = d.getNbCols();
2902 auto mat_shift = d.getMatShift();
2903
2904 // q0
2905 arrays_coord[0][4 * i + 0] = row;
2906 arrays_coord[1][4 * i + 0] = col;
2907 arrays_coord[2][4 * i + 0] = 0;
2908
2909 // q1
2910 arrays_coord[0][4 * i + 1] = row + nb_rows;
2911 arrays_coord[1][4 * i + 1] = col;
2912 arrays_coord[2][4 * i + 1] = 0;
2913
2914 // q2
2915 arrays_coord[0][4 * i + 2] = row + nb_rows;
2916 arrays_coord[1][4 * i + 2] = col + nb_cols;
2917 arrays_coord[2][4 * i + 2] = 0;
2918
2919 // q3
2920 arrays_coord[0][4 * i + 3] = row;
2921 arrays_coord[1][4 * i + 3] = col + nb_cols;
2922 arrays_coord[2][4 * i + 3] = 0;
2923
2924 // ele conn
2925 array[4 * i + 0] = starting_vertex + 4 * i + 0;
2926 array[4 * i + 1] = starting_vertex + 4 * i + 1;
2927 array[4 * i + 2] = starting_vertex + 4 * i + 2;
2928 array[4 * i + 3] = starting_vertex + 4 * i + 3;
2929
2930 auto prt = &(*block_mat_data->dataBlocksPtr)[d.getMatShift()];
2931 auto nrm2 = cblas_dnrm2(nb_rows * nb_cols, prt, 1);
2932 EntityHandle ele = starting_handle + i;
2933 CHKERR moab.tag_set_data(th_nrm2, &ele, 1, &nrm2);
2934 CHKERR moab.tag_set_data(th_mat_shift, &ele, 1, &mat_shift);
2935
2936 ++i;
2937 }
2938
2939 CHKERR iface->update_adjacencies(starting_handle, size, 4, array);
2940 CHKERR moab.write_file(filename.c_str(), "VTK", "");
2941
2943}
2944
2945} // namespace MoFEM
#define MOFEM_LOG_SEVERITY_SYNC(comm, severity)
Synchronise "SYNC" on curtain severity level.
#define MOFEM_TAG_AND_LOG(channel, severity, tag)
Tag and log in channel.
constexpr double a
@ COL
@ ROW
#define MAX_DOFS_ON_ENTITY
Maximal number of DOFs on entity.
#define CHK_THROW_MESSAGE(err, msg)
Check and throw MoFEM exception.
#define MoFEMFunctionReturnHot(a)
Last executable line of each PETSc function used for error handling. Replaces return()
@ NOSPACE
Definition definitions.h:83
#define MoFEMFunctionBegin
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
#define BITFIELDID_SIZE
max number of fields
@ MOFEM_OPERATION_UNSUCCESSFUL
Definition definitions.h:34
@ MOFEM_DATA_INCONSISTENCY
Definition definitions.h:31
@ MOFEM_NOT_IMPLEMENTED
Definition definitions.h:32
#define MoFEMFunctionReturn(a)
Last executable line of each PETSc function used for error handling. Replaces return()
#define CHKERR
Inline error check.
#define MoFEMFunctionBeginHot
First executable line of each MoFEM function, used for error handling. Final line of MoFEM functions ...
auto createDMVector(DM dm)
Get smart vector from DM.
Definition DMMoFEM.hpp:1234
PetscErrorCode DMoFEMLoopFiniteElementsUpAndLowRank(DM dm, const char fe_name[], MoFEM::FEMethod *method, int low_rank, int up_rank, CacheTupleWeakPtr cache_ptr=CacheTupleSharedPtr())
Executes FEMethod for finite elements in DM.
Definition DMMoFEM.cpp:557
#define MOFEM_LOG(channel, severity)
Log.
#define MOFEM_LOG_CHANNEL(channel)
Set and reset channel.
FTensor::Index< 'i', SPACE_DIM > i
const double c
speed of light (cm/ns)
const double v
phase velocity of light in medium (cm/ns)
const double n
refractive index of diffusive medium
static __CLPK_integer lapack_dgesv(__CLPK_integer n, __CLPK_integer nrhs, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_integer *ipiv, __CLPK_doublereal *b, __CLPK_integer ldb)
static __CLPK_integer lapack_dgetrf(__CLPK_integer m, __CLPK_integer n, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_integer *ipiv)
static __CLPK_integer lapack_dgetri(__CLPK_integer n, __CLPK_doublereal *a, __CLPK_integer lda, __CLPK_integer *ipiv, __CLPK_doublereal *work, __CLPK_integer lwork)
FTensor::Index< 'j', 3 > j
static MoFEMErrorCodeGeneric< PetscErrorCode > ierr
PetscErrorCode MoFEMErrorCode
MoFEM/PETSc error code.
uint128_t UId
Unique Id.
Definition Types.hpp:31
UBlasMatrix< double > MatrixDouble
Definition Types.hpp:77
UBlasVector< double > VectorDouble
Definition Types.hpp:68
UBlasVector< int > VectorInt
Definition Types.hpp:67
implementation of Data Operators for Forces and Sources
Definition Common.hpp:10
auto type_from_handle(const EntityHandle h)
get type from entity handle
MoFEMErrorCode MatSetValues< SchurElemMats >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2736
MoFEMErrorCode schur_mat_set_values_wrap(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2743
static auto cmp_uid_lo(const boost::weak_ptr< FieldEntity > &a, const UId &b)
static PetscErrorCode solve_add(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1325
OpSchurAssembleBase * createOpSchurAssembleEnd(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range > > field_ents, SmartPetscObj< AO > ao, SmartPetscObj< Mat > schur, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
Definition Schur.cpp:2663
static auto cmp_uid_hi(const UId &b, const boost::weak_ptr< FieldEntity > &a)
MoFEMErrorCode MatSetValues< SchurElemMatsBlock >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2812
auto isAllGather(IS is)
IS All gather.
static const bool debug
static MoFEMErrorCode mult_schur_block_shell(Mat mat, Vec x, Vec y, InsertMode iora, boost::function< int(DiagBlockIndex::BlockIndex::nth_index< 0 >::type::iterator)> shift_extractor, boost::shared_ptr< std::vector< double > > data_blocks_ptr, bool multiply_by_preconditioner, enum CBLAS_TRANSPOSE trans_A)
Definition Schur.cpp:1526
MoFEMErrorCode schurSwitchPreconditioner(boost::shared_ptr< BlockStructure > block_mat_data)
Switch preconditioner.
Definition Schur.cpp:2852
MoFEMErrorCode setSchurA00MatSolvePC(SmartPetscObj< PC > pc)
Set PC for A00 block.
Definition Schur.cpp:2705
constexpr bool debug_schur
Definition Schur.cpp:12
static PetscErrorCode solve_transpose_add(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1331
auto id_from_handle(const EntityHandle h)
MoFEMErrorCode MatSetValues< SchurElemMatsPreconditionedBlock >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2843
std::pair< SmartPetscObj< Mat >, boost::shared_ptr< BlockStructure > > SchurShellMatData
Definition Schur.hpp:98
std::pair< SmartPetscObj< Mat >, boost::shared_ptr< NestSchurData > > createSchurNestedMatrix(boost::shared_ptr< NestSchurData > schur_net_data_ptr)
Create a Mat Diag Blocks object.
Definition Schur.cpp:2630
static PetscErrorCode mult_transpose_add(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1309
MoFEMErrorCode MatSetValues< BlockStructure >(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2780
MoFEMErrorCode shell_block_mat_asmb_wrap(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2377
auto createISGeneral(MPI_Comm comm, PetscInt n, const PetscInt idx[], PetscCopyMode mode)
Creates a data structure for an index set containing a list of integers.
SmartPetscObj< Vec > vectorDuplicate(Vec vec)
Create duplicate vector of smart vector.
static PetscErrorCode solve(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1322
boost::shared_ptr< std::vector< double > > getBlockMatPrconditionerStorageMat(Mat B)
Get the Block Storage object.
Definition Schur.cpp:2217
static PetscErrorCode zero_rows_columns(Mat A, PetscInt N, const PetscInt rows[], PetscScalar diag, Vec x, Vec b)
Definition Schur.cpp:1335
boost::shared_ptr< BlockStructure > createBlockMatStructure(DM dm, SchurFEOpsFEandFields schur_fe_op_vec)
Create a Mat Diag Blocks object.
Definition Schur.cpp:1082
boost::shared_ptr< NestSchurData > createSchurNestedMatrixStruture(std::pair< SmartPetscObj< DM >, SmartPetscObj< DM > > dms, boost::shared_ptr< BlockStructure > block_mat_data_ptr, std::vector< std::string > fields_names, std::vector< boost::shared_ptr< Range > > field_ents, bool add_preconditioner_block)
Get the Schur Nest Mat Array object.
Definition Schur.cpp:2421
static PetscErrorCode mult_add(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1285
MoFEMErrorCode shell_block_mat_asmb_wrap_impl(BlockStructure *ctx, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora, boost::function< int(const DiagBlockIndex::Indexes *)> shift_extractor, boost::shared_ptr< std::vector< double > > data_blocks_ptr)
Definition Schur.cpp:2223
MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const double *ptr, InsertMode iora)
Assemble PETSc matrix.
static PetscErrorCode mult_transpose(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1297
auto field_bit_from_bit_number(const int bit_number)
get field bit id from bit number
OpSchurAssembleBase * createOpSchurZeroRowsAndCols(boost::shared_ptr< std::vector< unsigned char > > marker_ptr, double diag_val)
Create a Op Schur Zero Rows And Cols object.
Definition Schur.cpp:2675
static PetscErrorCode mat_zero(Mat m)
Definition Schur.cpp:1455
MoFEMErrorCode schurSaveBlockMesh(boost::shared_ptr< BlockStructure > block_mat_data, std::string filename)
Save block matrix as a mesh.
Definition Schur.cpp:2867
auto getInterfacePtr(DM dm)
Get the Interface Ptr object.
Definition DMMoFEM.hpp:1168
static PetscErrorCode mult(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1273
static PetscErrorCode solve_transpose(Mat mat, Vec x, Vec y)
Definition Schur.cpp:1328
static MoFEMErrorCode solve_schur_block_shell(Mat mat, Vec y, Vec x, InsertMode iora, enum CBLAS_TRANSPOSE trans_A)
Definition Schur.cpp:1734
std::vector< std::pair< std::string, std::vector< SchurFieldPair > > > SchurFEOpsFEandFields
Definition Schur.hpp:83
static MoFEMErrorCode setSchurBlockMatOps(Mat mat_raw)
Definition Schur.cpp:1467
MoFEMErrorCode shell_block_preconditioner_mat_asmb_wrap(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2397
auto getProblemPtr(DM dm)
get problem pointer from DM
Definition DMMoFEM.hpp:1179
SchurShellMatData createBlockMat(DM dm, boost::shared_ptr< BlockStructure > data)
Create a Schur Mat object.
Definition Schur.cpp:1492
MoFEMErrorCode assembleBlockMatSchur(MoFEM::Interface &m_field, Mat B, Mat S, std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range > > field_ents, SmartPetscObj< AO > ao)
Assemble Schur matrix.
Definition Schur.cpp:1895
OpSchurAssembleBase * createOpSchurAssembleBegin()
Definition Schur.cpp:2658
auto ent_form_type_and_id(const EntityType type, const EntityID id)
get entity handle from type and id
boost::shared_ptr< std::vector< double > > getBlockMatStorageMat(Mat B)
Get the Block Storage object.
Definition Schur.cpp:2210
constexpr AssemblyType A
SmartPetscObj< Mat > block_mat
FTensor::Index< 'm', 3 > m
const int N
Definition speed_test.cpp:3
bool multiplyByPreconditioner
Definition Schur.cpp:314
boost::shared_ptr< std::vector< double > > preconditionerBlocksPtr
Definition Schur.cpp:311
boost::shared_ptr< std::vector< double > > parentBlockStructurePtr
Definition Schur.cpp:312
boost::shared_ptr< std::vector< double > > dataBlocksPtr
Definition Schur.cpp:310
SmartPetscObj< Vec > ghostX
Definition Schur.cpp:307
SmartPetscObj< Vec > ghostY
Definition Schur.cpp:308
virtual FieldBitNumber get_field_bit_number(const std::string name) const =0
get field bit number
Deprecated interface functions.
block data indexes
Definition Schur.cpp:229
Indexes(UId uid_row, UId uid_col, UId uid_fe, int row, int col, int nb_rows, int nb_cols, int loc_row, int loc_col, int mat_shift)
Definition Schur.cpp:231
multi_index_container< Indexes, indexed_by< ordered_unique< composite_key< Indexes, const_mem_fun< Indexes, UId, &Indexes::getRowUId >, const_mem_fun< Indexes, UId, &Indexes::getColUId > > >, ordered_non_unique< const_mem_fun< Indexes, UId, &Indexes::getFEUId > >, ordered_non_unique< const_mem_fun< Indexes, UId, &Indexes::getRowUId > >, ordered_non_unique< const_mem_fun< Indexes, UId, &Indexes::getColUId > > > > BlockIndex
Definition Schur.cpp:300
virtual ~DiagBlockIndex()=default
BlockIndex blockIndex
blocks indexes storage
Definition Schur.cpp:302
static UId getLoFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
static UId getHiFieldEntityUId(const FieldBitNumber bit, const EntityHandle ent)
static UId getUniqueIdCalculate(const DofIdx dof, UId ent_uid)
Data on single entity (This is passed as argument to DataOperator::doWork)
const VectorFieldEntities & getFieldEntities() const
Get field entities (const version)
const VectorDofs & getFieldDofs() const
Get DOF data structures (const version)
const VectorInt & getIndices() const
Get global indices of degrees of freedom on entity.
[Storage and set boundary conditions]
Structure for user loop methods on finite elements.
static UId getHiBitNumberUId(const FieldBitNumber bit_number)
static auto getHandleFromUniqueId(const UId uid)
Get the Handle From Unique Id.
static auto getFieldBitNumberFromUniqueId(const UId uid)
Get the Field Bit Number From Unique Id.
UId getLocalUniqueIdCalculate()
Get the Local Unique Id Calculate object.
static UId getLoBitNumberUId(const FieldBitNumber bit_number)
static UId getLoLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
static UId getHiLocalEntityBitNumber(const char bit_number, const EntityHandle ent)
const FEMethod * getFEMethod() const
Return raw pointer to Finite Element Method object.
MoFEMErrorCode assembleSchurMat(Mat S, const UId &uid_row, VectorInt &row_ind, const UId &uid_col, VectorInt &col_ind, MatrixDouble &m, InsertMode iora)
Assemble Schur complement.
Definition Schur.cpp:38
Clear Schur complement internal data.
Definition Schur.cpp:61
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition Schur.cpp:523
Assemble Schur complement (Implementation)
Definition Schur.cpp:88
SmartPetscObj< AO > schurAO
Definition Schur.cpp:122
MoFEMErrorCode doWorkImpl(int side, EntityType type, EntitiesFieldData::EntData &data)
Definition Schur.cpp:594
SmartPetscObj< Mat > schurMat
Definition Schur.cpp:123
std::vector< boost::shared_ptr< Range > > fieldEnts
Definition Schur.cpp:121
OpSchurAssembleEndImpl(std::vector< std::string > fields_name, std::vector< boost::shared_ptr< Range > > field_ents, SmartPetscObj< AO > schur_ao, SmartPetscObj< Mat > schur_mat, bool sym_schur, bool symm_op)
Construct a new Op Schur Assemble End object.
Definition Schur.cpp:583
std::vector< std::string > fieldsName
Definition Schur.cpp:120
Assemble Schur complement.
Definition Schur.cpp:139
OpSchurZeroRowsAndCols(boost::shared_ptr< std::vector< unsigned char > > marker_ptr, double diag_val)
Definition Schur.cpp:536
boost::shared_ptr< std::vector< unsigned char > > markerPtr
Definition Schur.cpp:79
MoFEMErrorCode doWork(int side, EntityType type, EntitiesFieldData::EntData &data)
Operator for linear form, usually to calculate values on right hand side.
Definition Schur.cpp:542
static MatSetValuesPtr matSetValuesPreconditionedBlockPtr
Definition Schur.hpp:235
static MatSetValuesPtr matSetValuesPtr
backend assembly function
Definition Schur.hpp:231
static MatSetValuesPtr matSetValuesBlockPtr
backend assembly block mat function
Definition Schur.hpp:233
boost::function< MoFEMErrorCode(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)> MatSetValuesPtr
Definition Schur.hpp:230
static auto invertMat(MatrixDouble &m, MatrixDouble &inv)
Definition Schur.cpp:1039
static auto invertMat(MatrixDouble &m, MatrixDouble &inv)
Definition Schur.cpp:1003
static MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2799
static MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2831
Schur complement data storage.
Definition Schur.cpp:161
static MoFEMErrorCode assembleStorage(const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:347
static boost::ptr_vector< VectorInt > colIndices
Definition Schur.cpp:214
static boost::ptr_vector< MatrixDouble > locMats
Definition Schur.cpp:212
SchurElemMats(const size_t idx, const UId uid_row, const UId uid_col)
Definition Schur.cpp:342
virtual ~SchurElemMats()=default
static size_t maxIndexCounter
Definition Schur.cpp:216
static boost::ptr_vector< VectorInt > rowIndices
Definition Schur.cpp:213
auto & getMat() const
Definition Schur.cpp:170
static MoFEMErrorCode MatSetValues(Mat M, const EntitiesFieldData::EntData &row_data, const EntitiesFieldData::EntData &col_data, const MatrixDouble &mat, InsertMode iora)
Definition Schur.cpp:2757
auto & getColInd() const
Definition Schur.cpp:172
static boost::ptr_vector< SchurElemMats > schurElemMats
Definition Schur.cpp:215
auto & getRowInd() const
Definition Schur.cpp:171
static SchurElemStorage schurL2Storage
Definition Schur.cpp:218
const size_t iDX
Definition Schur.cpp:166
multi_index_container< const SchurElemMats *, indexed_by< ordered_unique< tag< uid_mi_tag >, composite_key< SchurElemMats, member< SchurElemMats, const UId, &SchurElemMats::uidRow >, member< SchurElemMats, const UId, &SchurElemMats::uidCol > > >, ordered_non_unique< tag< col_mi_tag >, member< SchurElemMats, const UId, &SchurElemMats::uidCol > > > > SchurElemStorage
Definition Schur.cpp:210
static PetscLogEvent MOFEM_EVENT_BlockStructureSetValues
Definition Schur.hpp:29
static PetscLogEvent MOFEM_EVENT_BlockStructureMult
Definition Schur.hpp:30
static PetscLogEvent MOFEM_EVENT_zeroRowsAndCols
Definition Schur.hpp:33
static PetscLogEvent MOFEM_EVENT_opSchurAssembleEnd
Definition Schur.hpp:28
static PetscLogEvent MOFEM_EVENT_BlockStructureSolve
Definition Schur.hpp:31
static PetscLogEvent MOFEM_EVENT_schurMatSetValues
Definition Schur.hpp:27
static PetscLogEvent MOFEM_EVENT_AssembleSchurMat
Definition Schur.hpp:32
intrusive_ptr for managing petsc objects