12 #include "../../config/config.hpp"
21 #include "moonolith_bounding_volume_with_span.hpp"
22 #include "moonolith_n_tree_mutator_factory.hpp"
23 #include "moonolith_n_tree_with_span_mutator_factory.hpp"
24 #include "moonolith_n_tree_with_tags_mutator_factory.hpp"
25 #include "moonolith_profiler.hpp"
26 #include "moonolith_redistribute.hpp"
27 #include "moonolith_sparse_matrix.hpp"
28 #include "moonolith_tree.hpp"
29 #include "par_moonolith.hpp"
33 using namespace mfem::internal;
38 struct ParMortarAssembler::Impl
42 std::shared_ptr<ParFiniteElementSpace>
source;
43 std::shared_ptr<ParFiniteElementSpace> destination;
44 std::vector<std::shared_ptr<MortarIntegrator>> integrators;
46 std::shared_ptr<HypreParMatrix> coupling_matrix;
47 std::shared_ptr<HypreParMatrix> mass_matrix;
49 std::shared_ptr<IterativeSolver> solver;
53 bool assemble_mass_and_coupling_together{
true};
58 solver = std::make_shared<BiCGSTABSolver>(comm);
59 solver->SetMaxIter(20000);
60 solver->SetRelTol(1e-6);
64 bool is_vector_fe()
const
66 bool is_vector_fe =
false;
67 for (
auto i_ptr : integrators)
69 if (i_ptr->is_vector_fe())
79 BilinearFormIntegrator * new_mass_integrator()
const
83 return new VectorFEMassIntegrator();
87 return new MassIntegrator();
96 impl_->solver = solver;
101 impl_->assemble_mass_and_coupling_together = value;
106 impl_->ensure_solver();
108 impl_->solver->SetMaxIter(max_solver_iterations);
112 const std::shared_ptr<MortarIntegrator> &integrator)
114 impl_->integrators.push_back(integrator);
119 impl_->verbose = verbose;
122 template <
int Dimension>
class ElementAdapter :
public moonolith::Serializable
125 using Bound = moonolith::AABBWithKDOPSpan<Dimension, double>;
126 using Point = moonolith::Vector<double, Dimension>;
128 inline int tag()
const {
return tag_; }
130 const Bound &bound()
const {
return bound_; }
132 Bound &bound() {
return bound_; }
134 void applyRW(moonolith::Stream &stream)
138 stream &element_handle_;
141 ElementAdapter(FiniteElementSpace &fe,
const long element,
142 const long element_handle,
const int tag)
143 : fe_(&fe), element_(element), element_handle_(element_handle), tag_(tag),
146 assert(element < fe.GetNE());
149 fe_->GetMesh()->GetPointMatrix(element, pts);
152 for (
int j = 0; j < pts.Width(); ++j)
154 for (
int i = 0; i < pts.Height(); ++i)
156 p[i] = pts.Elem(i, j);
159 bound_.static_bound() +=
p;
160 bound_.dynamic_bound() +=
p;
165 : fe_(nullptr), element_(-1), element_handle_(-1), tag_(-1),
168 inline long handle()
const {
return element_handle_; }
170 inline long element()
const {
return element_; }
172 inline const FiniteElement &
get()
const
175 assert(element_ < fe_->GetNE());
176 return *fe_->GetFE(element_);
179 inline const FiniteElementSpace &
space()
const
185 void set_dof_map(std::vector<long> *ptr) { dof_map_ = ptr; }
187 void get_elements_vdofs(Array<int> &vdofs)
const
189 fe_->GetElementVDofs(element_, vdofs);
193 assert(dof_map_->size() == vdofs.Size());
195 for (
int i = 0; i < vdofs.Size(); ++i)
197 vdofs[i] = dof_map_->at(i);
208 FiniteElementSpace *fe_;
210 long element_handle_;
213 std::vector<long> *dof_map_;
216 template <
int _Dimension>
class TreeTraits
219 enum { Dimension = _Dimension };
221 using Bound = moonolith::AABBWithKDOPSpan<Dimension, double>;
222 using DataType = mfem::ElementAdapter<Dimension>;
225 template <
int Dimension>
226 class MFEMTree :
public moonolith::Tree<TreeTraits<Dimension>>
229 using Traits = mfem::TreeTraits<Dimension>;
233 static std::shared_ptr<MFEMTree>
234 New(
const int maxElementsXNode = moonolith::DEFAULT_REFINE_MAX_ELEMENTS,
235 const int maxDepth = moonolith::DEFAULT_REFINE_DEPTH)
237 using namespace moonolith;
239 std::shared_ptr<MFEMTree> tree = std::make_shared<MFEMTree>();
240 std::shared_ptr<NTreeWithSpanMutatorFactory<MFEMTree>> factory =
241 std::make_shared<NTreeWithSpanMutatorFactory<MFEMTree>>();
242 factory->set_refine_params(maxElementsXNode, maxDepth);
243 tree->set_mutator_factory(factory);
247 static std::shared_ptr<MFEMTree>
248 New(
const std::shared_ptr<moonolith::Predicate> &predicate,
249 const int maxElementsXNode = moonolith::DEFAULT_REFINE_MAX_ELEMENTS,
250 const int maxDepth = moonolith::DEFAULT_REFINE_DEPTH)
252 using namespace moonolith;
256 return New(maxElementsXNode, maxDepth);
259 std::shared_ptr<MFEMTree> tree = std::make_shared<MFEMTree>();
260 std::shared_ptr<NTreeWithTagsMutatorFactory<MFEMTree>> factory =
261 std::make_shared<NTreeWithTagsMutatorFactory<MFEMTree>>(predicate);
262 factory->set_refine_params(maxElementsXNode, maxDepth);
263 tree->set_mutator_factory(factory);
268 class ElementDofMap :
public moonolith::Serializable
271 void read(moonolith::InputStream &is)
override
276 is.read(&global[0], n);
279 void write(moonolith::OutputStream &os)
const override
281 int n = global.size();
283 os.write(&global[0], n);
286 std::vector<long> global;
292 explicit Spaces(
const moonolith::Communicator &comm) : comm(comm)
294 must_destroy_attached[0] =
false;
295 must_destroy_attached[1] =
false;
298 Spaces(
const std::shared_ptr<ParFiniteElementSpace> &
source,
299 const std::shared_ptr<ParFiniteElementSpace> &destination)
302 spaces_.push_back(source);
303 spaces_.push_back(destination);
305 must_destroy_attached[0] =
false;
306 must_destroy_attached[1] =
false;
308 copy_global_dofs(*source, dof_maps_[0]);
309 copy_global_dofs(*destination, dof_maps_[1]);
315 FiniteElementCollection *fec =
nullptr;
317 for (
int i = 0; i < spaces_.size(); ++i)
319 if (spaces_[i] && must_destroy_attached[0])
321 m = spaces_[i]->GetMesh();
322 fec =
const_cast<FiniteElementCollection *
>(spaces_[i]->FEColl());
325 spaces_[i] = std::shared_ptr<FiniteElementSpace>();
333 inline long n_elements()
const
336 for (
auto s : spaces_)
347 inline std::vector<std::shared_ptr<FiniteElementSpace>> &spaces()
352 inline const std::vector<std::shared_ptr<FiniteElementSpace>> &
358 inline std::vector<ElementDofMap> &dof_map(
const int i)
365 inline const std::vector<ElementDofMap> &dof_map(
const int i)
const
372 inline void set_must_destroy_attached(
const int index,
const bool value)
376 must_destroy_attached[
index] = value;
380 std::vector<std::shared_ptr<FiniteElementSpace>> spaces_;
381 moonolith::Communicator comm;
382 std::vector<ElementDofMap> dof_maps_[2];
383 bool must_destroy_attached[2];
385 inline static void copy_global_dofs(ParFiniteElementSpace &fe,
386 std::vector<ElementDofMap> &dof_map)
388 dof_map.resize(fe.GetNE());
390 for (
int i = 0; i < fe.GetNE(); ++i)
392 fe.GetElementVDofs(i, vdofs);
393 for (
int k = 0; k < vdofs.Size(); ++k)
398 g_dof = fe.GetGlobalTDofNumber(vdofs[k]);
402 g_dof = -1 - fe.GetGlobalTDofNumber(-1 - vdofs[k]);
405 dof_map[i].global.push_back(g_dof);
411 template <
class Iterator>
412 static void write_space(
const Iterator &begin,
const Iterator &end,
413 FiniteElementSpace &
space,
414 const std::vector<ElementDofMap> &dof_map,
415 const int role, moonolith::OutputStream &os)
417 const int dim = space.GetMesh()->Dimension();
418 const long n_elements = std::distance(begin, end);
420 std::set<long> nodeIds;
421 std::map<long, long> mapping;
424 for (Iterator it = begin; it != end; ++it)
427 space.GetElementVertices(i, verts);
429 for (
int j = 0; j < verts.Size(); ++j)
431 nodeIds.insert(verts[j]);
435 long n_nodes = nodeIds.size();
438 os.request_space((n_elements * 8 + n_nodes * dim) *
439 (
sizeof(
double) +
sizeof(
long)));
441 auto fe_coll = space.FEColl();
442 const char *name = fe_coll->Name();
443 const int name_lenght = strlen(name);
447 os.write(name, name_lenght);
450 for (
auto nodeId : nodeIds)
452 mapping[nodeId] = index++;
461 for (
auto node_id : nodeIds)
463 double *v = space.GetMesh()->GetVertex(node_id);
464 for (
int i = 0; i <
dim; ++i)
471 for (Iterator it = begin; it != end; ++it)
474 space.GetElementVertices(k, verts);
476 const int attribute = space.GetAttribute(k);
477 const int e_n_nodes = verts.Size();
478 const int type = space.GetElementType(k);
479 const int order = space.GetOrder(k);
482 os << type << attribute << order << e_n_nodes;
484 for (
int i = 0; i < e_n_nodes; ++i)
486 auto it = mapping.find(verts[i]);
487 assert(it != mapping.end());
489 int index = it->second;
500 template <
class Iterator>
501 static void write_element_selection(
const Iterator &begin,
const Iterator &end,
502 const Spaces &spaces,
503 moonolith::OutputStream &os)
505 if (spaces.spaces().empty())
511 auto m = spaces.spaces()[0];
512 std::shared_ptr<FiniteElementSpace>
s =
nullptr;
514 if (spaces.spaces().size() > 1)
516 s = spaces.spaces()[1];
519 std::vector<long> source_selection;
520 std::vector<long> destination_selection;
522 bool met_destination_selection =
false;
524 for (Iterator it = begin; it != end; ++it)
528 if (m && index >= m->GetNE())
531 destination_selection.push_back(index);
535 met_destination_selection =
true;
536 destination_selection.push_back(index);
540 assert(!met_destination_selection);
541 assert(index < m->GetNE());
542 source_selection.push_back(index);
546 const bool has_source = !source_selection.empty();
547 const bool has_destination = !destination_selection.empty();
549 os << has_source << has_destination;
553 write_space(source_selection.begin(), source_selection.end(), *m,
554 spaces.dof_map(0), 0, os);
559 write_space(destination_selection.begin(), destination_selection.end(), *
s,
560 spaces.dof_map(1), 1, os);
564 static FiniteElementCollection *FECollFromName(
const std::string &comp_name)
569 static void read_space(moonolith::InputStream &is,
570 std::shared_ptr<FiniteElementSpace> &space,
571 std::vector<ElementDofMap> &dof_map)
576 int dim, role, name_lenght;
580 std::string name(name_lenght, 0);
581 is.read(&name[0], name_lenght);
591 auto fe_coll = FECollFromName(name);
592 auto mesh_ptr =
new Mesh(dim, n_nodes, n_elements);
594 for (
long i = 0; i < n_nodes; ++i)
597 for (
int i = 0; i <
dim; ++i)
603 mesh_ptr->AddVertex(v);
606 dof_map.resize(n_elements);
607 std::vector<int> e2v;
608 for (
long i = 0; i < n_elements; ++i)
611 int type, attribute, order, e_n_nodes;
612 is >> type >> attribute >> order >> e_n_nodes;
613 e2v.resize(e_n_nodes);
615 for (
int i = 0; i < e_n_nodes; ++i)
622 mesh_ptr->AddElement(NewElem(type, &e2v[0], attribute));
628 Finalize(*mesh_ptr,
true);
631 space = make_shared<FiniteElementSpace>(mesh_ptr, fe_coll);
634 static void read_spaces(moonolith::InputStream &is, Spaces &spaces)
636 bool has_source, has_destination;
637 is >> has_source >> has_destination;
639 spaces.spaces().resize(2);
643 read_space(is, spaces.spaces()[0], spaces.dof_map(0));
644 spaces.set_must_destroy_attached(0,
true);
648 spaces.spaces()[0] =
nullptr;
653 read_space(is, spaces.spaces()[1], spaces.dof_map(1));
654 spaces.set_must_destroy_attached(1,
true);
658 spaces.spaces()[1] =
nullptr;
662 template <
int Dimensions,
class Fun>
663 static bool Assemble(moonolith::Communicator &comm,
664 std::shared_ptr<ParFiniteElementSpace> &source,
665 std::shared_ptr<ParFiniteElementSpace> &destination,
666 Fun process_fun,
const moonolith::SearchSettings &settings,
669 using namespace moonolith;
671 typedef mfem::MFEMTree<Dimensions> NTreeT;
672 typedef typename NTreeT::DataContainer DataContainer;
673 typedef typename NTreeT::DataType Adapter;
675 long maxNElements = settings.max_elements;
676 long maxDepth = settings.max_depth;
678 const int n_elements_source = source->GetNE();
679 const int n_elements_destination = destination->GetNE();
680 const int n_elements = n_elements_source + n_elements_destination;
682 auto predicate = std::make_shared<MasterAndSlave>();
683 predicate->add(0, 1);
685 MOONOLITH_EVENT_BEGIN(
"create_adapters");
687 std::shared_ptr<NTreeT> tree = NTreeT::New(predicate, maxNElements, maxDepth);
688 tree->reserve(n_elements);
690 std::shared_ptr<Spaces> local_spaces =
691 std::make_shared<Spaces>(
source, destination);
696 for (
auto s : local_spaces->spaces())
700 for (
int i = 0; i < s->GetNE(); ++i)
702 Adapter
a(*s, i, offset + i, space_num);
703 a.set_dof_map(&local_spaces->dof_map(space_num)[i].global);
707 offset += s->GetNE();
713 tree->root()->bound().static_bound().enlarge(1e-6);
715 MOONOLITH_EVENT_END(
"create_adapters");
718 std::map<long, std::shared_ptr<Spaces>> spaces;
719 std::map<long, std::vector<std::shared_ptr<Spaces>>> migrated_spaces;
721 auto read = [&spaces, &migrated_spaces,
722 comm](
const long ownerrank,
const long ,
723 bool is_forwarding, DataContainer &data, InputStream &in)
725 CHECK_STREAM_READ_BEGIN(
"vol_proj", in);
727 std::shared_ptr<Spaces> proc_space = std::make_shared<Spaces>(comm);
729 read_spaces(in, *proc_space);
733 assert(!spaces[ownerrank]);
734 spaces[ownerrank] = proc_space;
738 migrated_spaces[ownerrank].push_back(proc_space);
741 data.reserve(data.size() + proc_space->n_elements());
745 for (
auto s : proc_space->spaces())
749 for (
int i = 0; i < s->GetNE(); ++i)
751 data.push_back(Adapter(*s, i, offset + i, space_num));
752 data.back().set_dof_map(&proc_space->dof_map(space_num)[i].global);
755 offset += s->GetNE();
761 CHECK_STREAM_READ_END(
"vol_proj", in);
764 auto write = [&local_spaces, &spaces,
765 &comm](
const long ownerrank,
const long ,
766 const std::vector<long>::const_iterator &begin,
767 const std::vector<long>::const_iterator &end,
768 const DataContainer &, OutputStream &
out)
770 CHECK_STREAM_WRITE_BEGIN(
"vol_proj",
out);
772 if (ownerrank == comm.rank())
774 write_element_selection(begin, end, *local_spaces,
out);
778 auto it = spaces.find(ownerrank);
779 assert(it != spaces.end());
780 std::shared_ptr<Spaces> spaceptr = it->second;
781 assert(std::distance(begin, end) > 0);
782 write_element_selection(begin, end, *spaceptr,
out);
785 CHECK_STREAM_WRITE_END(
"vol_proj",
out);
788 long n_false_positives = 0, n_intersections = 0;
789 auto fun = [&n_false_positives, &n_intersections,
790 &process_fun](Adapter &
source, Adapter &destination) ->
bool
792 bool ok = process_fun(source, destination);
808 moonolith::search_and_compute(comm, tree, predicate,
read,
write, fun,
813 long n_total_candidates = n_intersections + n_false_positives;
815 long n_collection[3] = {n_intersections, n_total_candidates,
818 comm.all_reduce(n_collection, 3, moonolith::MPISum());
822 mfem::out <<
"n_intersections: " << n_collection[0]
823 <<
", n_total_candidates: " << n_collection[1]
824 <<
", n_false_positives: " << n_collection[2] << std::endl;
831 static void add_matrix(
const Array<int> &destination_vdofs,
const Array<int> &source_vdofs,
832 const DenseMatrix &elem_mat, moonolith::SparseMatrix<double> &mat_buffer)
834 for (
int i = 0; i < destination_vdofs.Size(); ++i)
836 long dof_I = destination_vdofs[i];
846 for (
int j = 0; j < source_vdofs.Size(); ++j)
848 long dof_J = source_vdofs[j];
850 double sign_J = sign_I;
858 mat_buffer.add(dof_I, dof_J, sign_J * elem_mat.Elem(i, j));
864 const std::vector<moonolith::Integer> &destination_ranges,
867 moonolith::SparseMatrix<double> &mat_buffer) {
869 auto && comm = mat_buffer.comm();
871 moonolith::Redistribute<moonolith::SparseMatrix<double>> redist(comm.get());
872 redist.apply(destination_ranges, mat_buffer, moonolith::AddAssign<double>());
874 std::vector<int> I(s_offsets[1] - s_offsets[0] + 1);
877 std::vector<HYPRE_Int> J;
878 std::vector<double> data;
879 J.reserve(mat_buffer.n_local_entries());
880 data.reserve(J.size());
882 for (
auto it = mat_buffer.iter(); it; ++it)
884 I[it.row() - s_offsets[0] + 1]++;
885 J.push_back(it.col());
889 for (
int i = 1; i < I.size(); ++i)
894 return std::make_shared<HypreParMatrix>(
895 comm.get(), s_offsets[1] - s_offsets[0], mat_buffer.rows(),
896 mat_buffer.cols(), &I[0], &J[0], &data[0], s_offsets, m_offsets);
907 template <
int Dimensions>
909 Assemble(moonolith::Communicator &comm,
910 ParMortarAssembler::Impl &impl,
911 const moonolith::SearchSettings &settings,
const bool verbose)
913 auto &source = impl.source;
914 auto &destination = impl.destination;
916 auto &integrators = impl.integrators;
917 auto &pmat = impl.coupling_matrix;
920 bool lump_mass =
false;
923 for (
auto i_ptr : integrators)
925 max_q_order = std::max(i_ptr->GetQuadratureOrder(), max_q_order);
928 const int dim = source->GetMesh()->Dimension();
929 std::shared_ptr<Cut> cut =
NewCut(dim);
932 assert(
false &&
"NOT Supported!");
936 Array<int> source_vdofs, destination_vdofs;
938 DenseMatrix cumulative_elemmat;
939 IntegrationRule src_ir;
940 IntegrationRule dest_ir;
942 double local_element_matrices_sum = 0.0;
944 const auto m_global_n_dofs = source->GlobalVSize();
945 auto *m_offsets = source->GetDofOffsets();
947 const auto s_global_n_dofs = destination->GlobalVSize();
948 auto *s_offsets = destination->GetDofOffsets();
950 moonolith::SparseMatrix<double> mat_buffer(comm);
951 mat_buffer.set_size(s_global_n_dofs, m_global_n_dofs);
953 moonolith::SparseMatrix<double> mass_mat_buffer(comm);
954 mass_mat_buffer.set_size(s_global_n_dofs, s_global_n_dofs);
956 std::unique_ptr<BilinearFormIntegrator> mass_integr(impl.new_mass_integrator());
958 auto fun = [&](
const ElementAdapter<Dimensions> &
source,
959 const ElementAdapter<Dimensions> &destination) ->
bool
961 const auto &src = source.space();
962 const auto &dest = destination.space();
964 const int src_index = source.element();
965 const int dest_index = destination.element();
967 auto &src_fe = *src.GetFE(src_index);
968 auto &dest_fe = *dest.GetFE(dest_index);
970 ElementTransformation &dest_Trans =
971 *dest.GetElementTransformation(dest_index);
978 const int src_order = src_order_mult * src_fe.GetOrder();
979 const int dest_order = dest_order_mult * dest_fe.GetOrder();
981 int contraction_order = src_order + dest_order;
982 if(impl.assemble_mass_and_coupling_together) {
983 contraction_order = std::max(contraction_order, 2 * dest_order);
986 const int order = contraction_order + dest_order_mult * dest_Trans.OrderW() + max_q_order;
988 cut->SetIntegrationOrder(order);
990 if (cut->BuildQuadrature(src, src_index, dest, dest_index, src_ir,
994 ElementTransformation &src_Trans =
995 *src.GetElementTransformation(src_index);
997 source.get_elements_vdofs(source_vdofs);
998 destination.get_elements_vdofs(destination_vdofs);
1001 for (
auto i_ptr : integrators)
1005 i_ptr->AssembleElementMatrix(src_fe, src_ir, src_Trans, dest_fe,
1006 dest_ir, dest_Trans, cumulative_elemmat);
1011 i_ptr->AssembleElementMatrix(src_fe, src_ir, src_Trans, dest_fe,
1012 dest_ir, dest_Trans, elemmat);
1013 cumulative_elemmat += elemmat;
1017 local_element_matrices_sum += Sum(cumulative_elemmat);
1018 add_matrix(destination_vdofs, source_vdofs, cumulative_elemmat, mat_buffer);
1020 if(impl.assemble_mass_and_coupling_together) {
1021 mass_integr->SetIntRule(&dest_ir);
1022 mass_integr->AssembleElementMatrix(dest_fe, dest_Trans, elemmat);
1025 int n = destination_vdofs.Size();
1027 for(
int i = 0; i < n; ++i) {
1028 double row_sum = 0.;
1029 for(
int j = 0; j < n; ++j) {
1030 row_sum += elemmat(i, j);
1034 elemmat(i, i) = row_sum;
1038 add_matrix(destination_vdofs, destination_vdofs, elemmat, mass_mat_buffer);
1049 if (!Assemble<Dimensions>(comm, source, destination, fun, settings,
1057 double volumes[2] = {local_element_matrices_sum};
1058 comm.all_reduce(volumes, 2, moonolith::MPISum());
1064 mfem::out <<
"sum(B): " << volumes[0] << std::endl;
1068 std::vector<moonolith::Integer> destination_ranges(comm.size() + 1, 0);
1070 std::copy(s_offsets, s_offsets + 2, destination_ranges.begin() + comm.rank());
1071 comm.all_reduce(&destination_ranges[0], destination_ranges.size(),
1072 moonolith::MPIMax());
1080 if(impl.assemble_mass_and_coupling_together) {
1088 const std::shared_ptr<ParFiniteElementSpace> &source,
1089 const std::shared_ptr<ParFiniteElementSpace> &destination)
1092 impl_->comm = source->GetComm();
1094 impl_->destination = destination;
1099 assert(!impl_->integrators.empty() &&
1100 "it must have at least on integrator see class MortarIntegrator");
1102 moonolith::SearchSettings settings;
1108 moonolith::Communicator comm(impl_->comm);
1109 if (impl_->source->GetMesh()->Dimension() == 2)
1111 ok = mfem::Assemble<2>(comm, *impl_, settings,
1115 if (impl_->source->GetMesh()->Dimension() == 3)
1117 ok = mfem::Assemble<3>(comm, *impl_, settings,
1121 pmat = impl_->coupling_matrix;
1135 const bool verbose = impl_->verbose;
1137 if (!impl_->coupling_matrix)
1145 auto &B = *impl_->coupling_matrix;
1146 auto &D = *impl_->mass_matrix;
1148 auto &P_source = *impl_->source->Dof_TrueDof_Matrix();
1149 auto &P_destination = *impl_->destination->Dof_TrueDof_Matrix();
1151 impl_->ensure_solver();
1156 impl_->solver->SetOperator(D);
1159 impl_->solver->SetPrintLevel(3);
1162 Vector P_x_src_fun(B.Width());
1163 P_source.MultTranspose(src_fun, P_x_src_fun);
1165 Vector B_x_src_fun(B.Height());
1168 B.Mult(P_x_src_fun, B_x_src_fun);
1170 Vector R_x_dest_fun(D.Height());
1173 impl_->solver->Mult(B_x_src_fun, R_x_dest_fun);
1175 P_destination.Mult(R_x_dest_fun, dest_fun);
1180 int converged = impl_->solver->GetFinalNorm() < 1e-5;
1181 MPI_Allreduce(MPI_IN_PLACE, &converged, 1, MPI_INT, MPI_MIN, impl_->comm);
1187 using namespace std;
1188 const bool verbose = impl_->verbose;
1189 static const bool dof_transformation =
true;
1191 moonolith::Communicator comm(impl_->comm);
1196 moonolith::root_describe(
1197 "--------------------------------------------------------"
1198 "\nAssembly begin: ",
1202 if (!
Assemble(impl_->coupling_matrix))
1209 moonolith::root_describe(
1211 "--------------------------------------------------------",
1215 if (dof_transformation)
1217 impl_->coupling_matrix.reset(
RAP(impl_->destination->Dof_TrueDof_Matrix(),
1218 impl_->coupling_matrix.get(),
1219 impl_->source->Dof_TrueDof_Matrix()));
1222 auto &B = *impl_->coupling_matrix;
1224 if(!impl_->assemble_mass_and_coupling_together) {
1229 impl_->mass_matrix = shared_ptr<HypreParMatrix>(b_form.ParallelAssemble());
1235 if (verbose && comm.is_root())
1237 mfem::out <<
"P in R^(" << impl_->source->Dof_TrueDof_Matrix()->Height();
1238 mfem::out <<
" x " << impl_->source->Dof_TrueDof_Matrix()->Width() <<
")\n";
1240 << impl_->destination->Dof_TrueDof_Matrix()->Height();
1241 mfem::out <<
" x " << impl_->destination->Dof_TrueDof_Matrix()->Width()
1245 if (dof_transformation && impl_->assemble_mass_and_coupling_together)
1247 impl_->mass_matrix.reset(
RAP(impl_->destination->Dof_TrueDof_Matrix(),
1248 impl_->mass_matrix.get(),
1249 impl_->destination->Dof_TrueDof_Matrix()));
1252 auto &D = *impl_->mass_matrix;
1256 if (verbose && comm.is_root())
1258 mfem::out <<
"--------------------------------------------------------"
1260 mfem::out <<
"B in R^(" << B.GetGlobalNumRows() <<
" x "
1261 << B.GetGlobalNumCols() <<
")" << std::endl;
1262 mfem::out <<
"D in R^(" << D.GetGlobalNumRows() <<
" x "
1263 << D.GetGlobalNumCols() <<
")" << std::endl;
1264 mfem::out <<
"--------------------------------------------------------"
1278 double sum_Dv = Dv.
Sum();
1279 comm.all_reduce(&sum_Dv, 1, MPI_SUM);
1283 mfem::out <<
"sum(D): " << sum_Dv << std::endl;
1292 #endif // MFEM_USE_MPI
void SetVerbose(const bool verbose)
Expose process details with verbose output.
int order_multiplier(const Geometry::Type type, const int dim)
bool Apply(const ParGridFunction &src_fun, ParGridFunction &dest_fun)
transfer a function from source to destination. It requires that the Update function is called before...
bool Transfer(const ParGridFunction &src_fun, ParGridFunction &dest_fun)
transfer a function from source to destination. if the transfer is to be performed multiple times use...
void AddMortarIntegrator(const std::shared_ptr< MortarIntegrator > &integrator)
This method must be called before Assemble or Transfer. It will assemble the operator in all intersec...
void source(const Vector &x, Vector &f)
void write(std::ostream &os, T value)
Write 'value' to stream.
virtual void Update()
Transform by the Space UpdateMatrix (e.g., on Mesh change).
ParMortarAssembler(const std::shared_ptr< ParFiniteElementSpace > &source, const std::shared_ptr< ParFiniteElementSpace > &destination)
constructs the object with source and destination spaces
T read(std::istream &is)
Read a value from the stream and return it.
void RAP(const DenseMatrix &A, const DenseMatrix &P, DenseMatrix &RAP)
std::shared_ptr< HypreParMatrix > convert_to_hypre_matrix(const std::vector< moonolith::Integer > &destination_ranges, HYPRE_BigInt *s_offsets, HYPRE_BigInt *m_offsets, moonolith::SparseMatrix< double > &mat_buffer)
double p(const Vector &x, double t)
std::shared_ptr< Cut > NewCut(const int dim)
void SetMaxSolverIterations(const int max_solver_iterations)
Control the maximum numbers of conjugate gradients steps for mass matrix inversion.
static FiniteElementCollection * New(const char *name)
Factory method: return a newly allocated FiniteElementCollection according to the given name...
bool Assemble(std::shared_ptr< HypreParMatrix > &B)
assembles the coupling matrix B. B : source -> destination If u is a coefficient associated with sour...
int index(int i, int j, int nx, int ny)
void SetAssembleMassAndCouplingTogether(const bool value)
Control if the Mass matrix is computed together with the coupling operator every time.
void pts(int iphi, int t, double x[])
void SetSolver(const std::shared_ptr< IterativeSolver > &solver)
Changes the solver to be used for solving the mass-matrix linear system.
Class for parallel grid function.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
bool Update()
assembles the various components necessary for the transfer. To be called before calling the Apply fu...
double Sum() const
Return the sum of the vector entries.