preCICE v3.2.0
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
ParticipantImpl.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <Eigen/src/Core/util/Meta.h>
3#include <algorithm>
4#include <array>
5#include <cmath>
6#include <deque>
7#include <filesystem>
8#include <functional>
9#include <iterator>
10#include <memory>
11#include <numeric>
12#include <optional>
13#include <ostream>
14#include <sstream>
15#include <tuple>
16#include <utility>
17
18#include "ParticipantImpl.hpp"
20#include "com/Communication.hpp"
21#include "com/SharedPointer.hpp"
25#include "io/Export.hpp"
26#include "io/ExportContext.hpp"
27#include "io/SharedPointer.hpp"
29#include "logging/LogMacros.hpp"
30#include "m2n/BoundM2N.hpp"
31#include "m2n/M2N.hpp"
32#include "m2n/SharedPointer.hpp"
34#include "mapping/Mapping.hpp"
38#include "math/differences.hpp"
39#include "math/geometry.hpp"
40#include "mesh/Data.hpp"
41#include "mesh/Edge.hpp"
42#include "mesh/Mesh.hpp"
44#include "mesh/Utils.hpp"
45#include "mesh/Vertex.hpp"
64#include "precice/impl/versions.hpp"
65#include "profiling/Event.hpp"
69#include "utils/EigenIO.hpp"
70#include "utils/Helpers.hpp"
71#include "utils/IntraComm.hpp"
72#include "utils/Parallel.hpp"
73#include "utils/Petsc.hpp"
74#include "utils/algorithm.hpp"
75#include "utils/assertion.hpp"
76#include "xml/XMLTag.hpp"
77
79
80namespace precice::impl {
81
83 std::string_view participantName,
84 std::string_view configurationFileName,
85 int solverProcessIndex,
86 int solverProcessSize,
87 std::optional<void *> communicator)
88 : _accessorName(participantName),
89 _accessorProcessRank(solverProcessIndex),
90 _accessorCommunicatorSize(solverProcessSize)
91{
92
94 "This participant's name is an empty string. "
95 "When constructing a preCICE interface you need to pass the name of the "
96 "participant as first argument to the constructor.");
98
99 PRECICE_CHECK(!communicator || communicator.value() != nullptr,
100 "Passing \"nullptr\" as \"communicator\" to Participant constructor is not allowed. "
101 "Please use the Participant constructor without the \"communicator\" argument, if you don't want to pass an MPI communicator.");
103 "The solver process index needs to be a non-negative number, not: {}. "
104 "Please check the value given when constructing a preCICE interface.",
107 "The solver process size needs to be a positive number, not: {}. "
108 "Please check the value given when constructing a preCICE interface.",
111 "The solver process index, currently: {} needs to be smaller than the solver process size, currently: {}. "
112 "Please check the values given when constructing a preCICE interface.",
114
117 Event e("construction", profiling::Fundamental);
118
119 // Set the global communicator to the passed communicator.
120 // This is a noop if preCICE is not configured with MPI.
121#ifndef PRECICE_NO_MPI
122 Event e3("com.initializeMPI", profiling::Fundamental);
123 if (communicator.has_value()) {
124 auto commptr = static_cast<utils::Parallel::Communicator *>(communicator.value());
126 } else {
128 }
129
130 {
131 const auto currentRank = utils::Parallel::current()->rank();
132 PRECICE_CHECK(_accessorProcessRank == currentRank,
133 "The solver process index given in the preCICE interface constructor({}) does not match the rank of the passed MPI communicator ({}).",
134 _accessorProcessRank, currentRank);
135 const auto currentSize = utils::Parallel::current()->size();
137 "The solver process size given in the preCICE interface constructor({}) does not match the size of the passed MPI communicator ({}).",
138 _accessorCommunicatorSize, currentSize);
139 }
140 e3.stop();
141#endif
142
143 Event e1("configure", profiling::Fundamental);
144 configure(configurationFileName);
145 e1.stop();
146
147 // Backend settings have been configured
148 Event e2("startProfilingBackend");
150 e2.stop();
151
152 PRECICE_DEBUG("Initialize intra-participant communication");
155 }
156
157 e.stop();
159}
160
162{
163 if (_state != State::Finalized) {
164 PRECICE_INFO("Implicitly finalizing in destructor");
165 finalize();
166 }
167}
168
170 std::string_view configurationFileName)
171{
172
179 _configHash = xml::configure(config.getXMLTag(), context, configurationFileName);
180 if (_accessorProcessRank == 0) {
181 PRECICE_INFO("This is preCICE version {}", PRECICE_VERSION);
182 PRECICE_INFO("Revision info: {}", precice::preciceRevision);
183 PRECICE_INFO("Build type: "
184#ifndef NDEBUG
185 "Debug"
186#else // NDEBUG
187 "Release"
188#ifndef PRECICE_NO_DEBUG_LOG
189 " + debug log"
190#else
191 " (without debug log)"
192#endif
193#ifndef PRECICE_NO_TRACE_LOG
194 " + trace log"
195#endif
196#ifndef PRECICE_NO_ASSERTIONS
197 " + assertions"
198#endif
199#endif // NDEBUG
200 );
201 try {
202 PRECICE_INFO("Working directory \"{}\"", std::filesystem::current_path().string());
203 } catch (std::filesystem::filesystem_error &fse) {
204 PRECICE_INFO("Working directory unknown due to error \"{}\"", fse.what());
205 }
206 PRECICE_INFO("Configuring preCICE with configuration \"{}\"", configurationFileName);
207 PRECICE_INFO("I am participant \"{}\"", _accessorName);
208 }
209
211
213
218 _accessor->setMeshIdManager(config.getMeshConfiguration()->extractMeshIdManager());
219
220 PRECICE_ASSERT(_accessorCommunicatorSize == 1 || _accessor->useIntraComm(),
221 "A parallel participant needs an intra-participant communication");
222 PRECICE_CHECK(not(_accessorCommunicatorSize == 1 && _accessor->useIntraComm()),
223 "You cannot use an intra-participant communication with a serial participant. "
224 "If you do not know exactly what an intra-participant communication is and why you want to use it "
225 "you probably just want to remove the intraComm tag from the preCICE configuration.");
226
228
229 _participants = config.getParticipantConfiguration()->getParticipants();
231
232 PRECICE_CHECK(_participants.size() > 1,
233 "In the preCICE configuration, only one participant is defined. "
234 "One participant makes no coupled simulation. "
235 "Please add at least another one.");
237
240 _couplingScheme = cplSchemeConfig->getCouplingScheme(_accessorName);
241
242 // Register all MeshIds to the lock, but unlock them straight away as
243 // writing is allowed after configuration.
244 for (const MeshContext *meshContext : _accessor->usedMeshContexts()) {
245 _meshLock.add(meshContext->mesh->getName(), false);
246 }
247}
248
250{
252 PRECICE_CHECK(_state != State::Finalized, "initialize() cannot be called after finalize().");
253 PRECICE_CHECK(_state != State::Initialized, "initialize() may only be called once.");
254 PRECICE_ASSERT(not _couplingScheme->isInitialized());
255
257 PRECICE_CHECK(not failedToInitialize,
258 "Initial data has to be written to preCICE before calling initialize(). "
259 "After defining your mesh, call requiresInitialData() to check if the participant is required to write initial data using the writeData() function.");
260
261 // Enforce that all user-created events are stopped to prevent incorrect nesting.
262 PRECICE_CHECK(_userEvents.empty(), "There are unstopped user defined events. Please stop them using stopLastProfilingSection() before calling initialize().");
263
264 _solverInitEvent.reset();
266
267 for (const auto &context : _accessor->usedMeshContexts()) {
268 if (context->provideMesh) {
269 e.addData("meshSize" + context->mesh->getName(), context->mesh->nVertices());
270 }
271 }
272
274 setupWatcher();
275
277
278 for (auto &context : _accessor->writeDataContexts()) {
279 const double startTime = 0.0;
280 context.storeBufferedData(startTime);
281 }
282
285
286 PRECICE_DEBUG("Initialize coupling schemes");
287 Event e1("initalizeCouplingScheme", profiling::Fundamental);
288 _couplingScheme->initialize();
289 e1.stop();
290
293
295
297
298 e.stop();
299
301 PRECICE_INFO(_couplingScheme->printCouplingState());
303}
304
306{
309 PRECICE_INFO("Reinitializing Participant");
310 Event e("reinitialize", profiling::Fundamental);
312
313 for (const auto &context : _accessor->usedMeshContexts()) {
314 if (context->provideMesh) {
315 e.addData("meshSize" + context->mesh->getName(), context->mesh->nVertices());
316 }
317 }
318
320 setupWatcher();
321
322 PRECICE_DEBUG("Reinitialize coupling schemes");
323 _couplingScheme->reinitialize();
324}
325
327{
329
330 // TODO only preprocess changed meshes
331 PRECICE_DEBUG("Preprocessing provided meshes");
332 for (MeshContext *meshContext : _accessor->usedMeshContexts()) {
333 if (meshContext->provideMesh) {
334 auto &mesh = *(meshContext->mesh);
335 Event e("preprocess." + mesh.getName());
336 meshContext->mesh->preprocess();
337 }
338 }
339
340 // Setup communication
341
342 PRECICE_INFO("Setting up primary communication to coupling partner/s");
343 Event e2("connectPrimaries");
344 for (auto &m2nPair : _m2ns) {
345 auto &bm2n = m2nPair.second;
346 bool requesting = bm2n.isRequesting;
347 if (bm2n.m2n->isConnected()) {
348 PRECICE_DEBUG("Primary connection {} {} already connected.", (requesting ? "from" : "to"), bm2n.remoteName);
349 } else {
350 PRECICE_DEBUG((requesting ? "Awaiting primary connection from {}" : "Establishing primary connection to {}"), bm2n.remoteName);
351 bm2n.prepareEstablishment();
352 bm2n.connectPrimaryRanks(_configHash);
353 PRECICE_DEBUG("Established primary connection {} {}", (requesting ? "from " : "to "), bm2n.remoteName);
354 }
355 }
356 e2.stop();
357
358 PRECICE_INFO("Primary ranks are connected");
359
360 Event e3("repartitioning");
361 // clears the mappings as well (see clearMappings)
363
364 PRECICE_INFO("Setting up preliminary secondary communication to coupling partner/s");
365 for (auto &m2nPair : _m2ns) {
366 auto &bm2n = m2nPair.second;
367 bm2n.preConnectSecondaryRanks();
368 }
369
371 e3.stop();
372
373 PRECICE_INFO("Setting up secondary communication to coupling partner/s");
374 Event e4("connectSecondaries");
375 for (auto &m2nPair : _m2ns) {
376 auto &bm2n = m2nPair.second;
377 bm2n.connectSecondaryRanks();
378 PRECICE_DEBUG("Established secondary connection {} {}", (bm2n.isRequesting ? "from " : "to "), bm2n.remoteName);
379 }
380 PRECICE_INFO("Secondary ranks are connected");
381
382 for (auto &m2nPair : _m2ns) {
383 m2nPair.second.cleanupEstablishment();
384 }
385}
386
388{
390 PRECICE_DEBUG("Initialize watchpoints");
391 for (PtrWatchPoint &watchPoint : _accessor->watchPoints()) {
392 watchPoint->initialize();
393 }
394 for (PtrWatchIntegral &watchIntegral : _accessor->watchIntegrals()) {
395 watchIntegral->initialize();
396 }
397}
398
400 double computedTimeStepSize)
401{
402
403 PRECICE_TRACE(computedTimeStepSize);
404
405 // Enforce that all user-created events are stopped to prevent incorrect nesting.
406 PRECICE_CHECK(_userEvents.empty(), "There are unstopped user defined events. Please stop them using stopLastProfilingSection() before calling advance().");
407
408 // Events for the solver time, stopped when we enter, restarted when we leave advance
409 PRECICE_ASSERT(_solverAdvanceEvent, "The advance event is created in initialize");
410 _solverAdvanceEvent->stop();
411
413
414 PRECICE_CHECK(_state != State::Constructed, "initialize() has to be called before advance().");
415 PRECICE_CHECK(_state != State::Finalized, "advance() cannot be called after finalize().");
416 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before advance().");
417 PRECICE_ASSERT(_couplingScheme->isInitialized());
418 PRECICE_CHECK(isCouplingOngoing(), "advance() cannot be called when isCouplingOngoing() returns false.");
419 PRECICE_CHECK(!math::equals(computedTimeStepSize, 0.0), "advance() cannot be called with a time step size of 0.");
420 PRECICE_CHECK(computedTimeStepSize > 0.0, "advance() cannot be called with a negative time step size {}.", computedTimeStepSize);
422
423#ifndef NDEBUG
424 PRECICE_DEBUG("Synchronize time step size");
426 syncTimestep(computedTimeStepSize);
427 }
428#endif
429
430 // Update the coupling scheme time state. Necessary to get correct remainder.
431 const bool isAtWindowEnd = _couplingScheme->addComputedTime(computedTimeStepSize);
432
433 if (_allowsRemeshing) {
434 if (isAtWindowEnd) {
435 auto totalMeshChanges = getTotalMeshChanges();
436 clearStamplesOfChangedMeshes(totalMeshChanges);
437
438 int sumOfChanges = std::accumulate(totalMeshChanges.begin(), totalMeshChanges.end(), 0);
439 if (reinitHandshake(sumOfChanges)) {
440 reinitialize();
441 }
442 } else {
443 PRECICE_CHECK(_meshLock.checkAll(), "The time window needs to end after remeshing.");
444 }
445 }
446
447 const double timeSteppedTo = _couplingScheme->getTime();
448 const auto dataToReceive = _couplingScheme->implicitDataToReceive();
449
450 handleDataBeforeAdvance(isAtWindowEnd, timeSteppedTo);
451
453
454 // In clase if an implicit scheme, this may be before timeSteppedTo
455 const double timeAfterAdvance = _couplingScheme->getTime();
456 const bool timeWindowComplete = _couplingScheme->isTimeWindowComplete();
457
458 handleDataAfterAdvance(isAtWindowEnd, timeWindowComplete, timeSteppedTo, timeAfterAdvance, dataToReceive);
459
460 PRECICE_INFO(_couplingScheme->printCouplingState());
461
462 PRECICE_DEBUG("Mapped {} samples in write mappings and {} samples in read mappings",
464
466
467 e.stop();
468 _solverAdvanceEvent->start();
469}
470
471void ParticipantImpl::handleDataBeforeAdvance(bool reachedTimeWindowEnd, double timeSteppedTo)
472{
473 // We only have to care about write data, in case substeps are enabled
474 // OR we are at the end of a timewindow, otherwise, we simply erase
475 // them as they have no relevance for the coupling (without time
476 // interpolation, only the time window end is relevant), the resetting
477 // happens regardless of the if-condition.
478 if (reachedTimeWindowEnd || _couplingScheme->requiresSubsteps()) {
479
480 // Here, we add the written data to the waveform storage. In the
481 // mapWrittenData, we then take samples from the storage and execute
482 // the mapping using waveform samples on the (for write mappings) "to"
483 // side.
484 samplizeWriteData(timeSteppedTo);
485 }
486
488
489 // Reset mapping counters here to cover subcycling
492
493 if (reachedTimeWindowEnd) {
494 mapWrittenData(_couplingScheme->getTimeWindowStart());
496 }
497}
498
499void ParticipantImpl::handleDataAfterAdvance(bool reachedTimeWindowEnd, bool isTimeWindowComplete, double timeSteppedTo, double timeAfterAdvance, const cplscheme::ImplicitData &receivedData)
500{
501 if (!reachedTimeWindowEnd) {
502 // We are subcycling
503 return;
504 }
505
507 // Move to next time window
508 PRECICE_ASSERT(math::greaterEquals(timeAfterAdvance, timeSteppedTo), "We must have stayed or moved forwards in time (min-time-step-size).", timeAfterAdvance, timeSteppedTo);
509
510 // As we move forward, there may now be old samples lying around
511 trimOldDataBefore(_couplingScheme->getTimeWindowStart());
512 } else {
513 // We are iterating
514 PRECICE_ASSERT(math::greater(timeSteppedTo, timeAfterAdvance), "We must have moved back in time!");
515
516 trimSendDataAfter(timeAfterAdvance);
517 }
518
519 if (reachedTimeWindowEnd) {
520 trimReadMappedData(timeAfterAdvance, isTimeWindowComplete, receivedData);
521 mapReadData();
523 }
524
525 // Required for implicit coupling
526 for (auto &context : _accessor->readDataContexts()) {
527 context.invalidateMappingCache();
528 }
529
530 // Strictly speaking, the write direction is not relevant here, but we will add it for the sake of completenss
531 for (auto &context : _accessor->writeDataContexts()) {
532 context.invalidateMappingCache();
533 }
534
536 // Reset initial guesses for iterative mappings
537 for (auto &context : _accessor->readDataContexts()) {
538 context.resetInitialGuesses();
539 }
540 for (auto &context : _accessor->writeDataContexts()) {
541 context.resetInitialGuesses();
542 }
543 }
544
546}
547
549{
550 // store buffered write data in sample storage and reset the buffer
551 for (auto &context : _accessor->writeDataContexts()) {
552
553 // Finalize conservative write mapping, later we reset
554 // the buffer in resetWrittenData
555
556 // Note that "samplizeWriteData" operates on _providedData of the
557 // DataContext, which is for just-in-time mappings the data we write
558 // on the received mesh.
559 // For just-in-time mappings, the _providedData should contain by now
560 // the "just-in-time" mapped data. However, it would be wasteful to
561 // execute expensive parts (in particular solving the RBF systems)
562 // for each writeAndMapData call. Thus, we create a DataCache during
563 // the writeAndMapData API calls, which contains pre-processed data
564 // values. Here, we now need to finalize the just-in-time mappings,
565 // before we can add it to the waveform buffer.
566 // For now, this only applies to just-in-time write mappings
567
568 context.completeJustInTimeMapping();
569 context.storeBufferedData(time);
570 }
571}
572
574{
575 for (auto &context : _accessor->usedMeshContexts()) {
576 for (const auto &name : context->mesh->availableData()) {
577 context->mesh->data(name)->timeStepsStorage().trimBefore(time);
578 }
579 }
580}
581
583{
584 for (auto &context : _accessor->writeDataContexts()) {
585 context.trimAfter(time);
586 }
587}
588
590{
592 PRECICE_CHECK(_state != State::Finalized, "finalize() may only be called once.");
593
594 // First we gracefully stop all existing user events and finally the last solver.advance event
595 while (!_userEvents.empty()) {
596 // Ensure reverse destruction order for correct nesting
597 _userEvents.pop_back();
598 }
599 _solverAdvanceEvent.reset();
600
601 Event e("finalize", profiling::Fundamental);
602
603 if (_state == State::Initialized) {
604
605 PRECICE_ASSERT(_couplingScheme->isInitialized());
606 PRECICE_DEBUG("Finalize coupling scheme");
607 _couplingScheme->finalize();
608
610 }
611
612 // Release ownership
614 _participants.clear();
616
617 // Close Connections
618 PRECICE_DEBUG("Close intra-participant communication");
620 utils::IntraComm::getCommunication()->closeConnection();
622 }
623 _m2ns.clear();
624
625 // Stop and print Event logging
626 e.stop();
627
628 // Finalize PETSc and Events first
630// This will lead to issues if we call finalize afterwards again
631#ifndef PRECICE_NO_GINKGO
633#endif
635
636 // Finally clear events and finalize MPI
639}
640
642{
643 PRECICE_TRACE(meshName);
645 return _accessor->usedMeshContext(meshName).mesh->getDimensions();
646}
647
649{
650 PRECICE_TRACE(meshName, dataName);
652 PRECICE_VALIDATE_DATA_NAME(meshName, dataName);
653 return _accessor->usedMeshContext(meshName).mesh->data(dataName)->getDimensions();
654}
655
657{
659 PRECICE_CHECK(_state != State::Finalized, "isCouplingOngoing() cannot be called after finalize().");
660 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before isCouplingOngoing() can be evaluated.");
661 return _couplingScheme->isCouplingOngoing();
662}
663
665{
667 PRECICE_CHECK(_state != State::Constructed, "initialize() has to be called before isTimeWindowComplete().");
668 PRECICE_CHECK(_state != State::Finalized, "isTimeWindowComplete() cannot be called after finalize().");
669 return _couplingScheme->isTimeWindowComplete();
670}
671
673{
674 PRECICE_CHECK(_state != State::Finalized, "getMaxTimeStepSize() cannot be called after finalize().");
675 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before getMaxTimeStepSize() can be evaluated.");
676 const double nextTimeStepSize = _couplingScheme->getNextTimeStepMaxSize();
677 // PRECICE_ASSERT(!math::equals(nextTimeStepSize, 0.0), nextTimeStepSize); // @todo requires https://github.com/precice/precice/issues/1904
678 // PRECICE_ASSERT(math::greater(nextTimeStepSize, 0.0), nextTimeStepSize); // @todo requires https://github.com/precice/precice/issues/1904
679
680 // safeguard needed because _couplingScheme->getNextTimeStepMaxSize() returns 0, if not isCouplingOngoing()
681 // actual case where we want to warn the user
683 isCouplingOngoing() && not math::greater(nextTimeStepSize, 0.0, 100 * math::NUMERICAL_ZERO_DIFFERENCE),
684 "preCICE just returned a maximum time step size of {}. Such a small value can happen if you use many substeps per time window over multiple time windows due to added-up differences of machine precision.",
685 nextTimeStepSize);
686 return nextTimeStepSize;
687}
688
690{
692 PRECICE_CHECK(_state == State::Constructed, "requiresInitialData() has to be called before initialize().");
694 if (required) {
696 }
697 return required;
698}
699
701{
703 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before requiresWritingCheckpoint().");
705 if (required) {
707 }
708 return required;
709}
710
712{
714 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before requiresReadingCheckpoint().");
716 if (required) {
718 }
719 return required;
720}
721
723{
725 MeshContext &context = _accessor->usedMeshContext(meshName);
727}
728
730 std::string_view dataName) const
731{
732 PRECICE_VALIDATE_DATA_NAME(meshName, dataName);
733 // Read data never requires gradients
734 if (!_accessor->isDataWrite(meshName, dataName))
735 return false;
736
737 WriteDataContext &context = _accessor->writeDataContext(meshName, dataName);
738 return context.hasGradient();
739}
740
742 std::string_view meshName) const
743{
744 PRECICE_TRACE(meshName);
745 PRECICE_REQUIRE_MESH_USE(meshName);
746 // In case we access received mesh data: check, if the requested mesh data has already been received.
747 // Otherwise, the function call doesn't make any sense
748 PRECICE_CHECK((_state == State::Initialized) || _accessor->isMeshProvided(meshName), "initialize() has to be called before accessing"
749 " data of the received mesh \"{}\" on participant \"{}\".",
750 meshName, _accessor->getName());
751 MeshContext &context = _accessor->usedMeshContext(meshName);
752 PRECICE_ASSERT(context.mesh.get() != nullptr);
753
754 // Returns true if we have api access configured and we run in parallel and have a received mesh
755 if ((context.userDefinedAccessRegion || requiresUserDefinedAccessRegion(meshName)) && _accessor->isDirectAccessAllowed(meshName)) {
756 // filter nVertices to the actual number of vertices queried by the user
757 PRECICE_CHECK(context.userDefinedAccessRegion, "The function getMeshVertexSize was called on the received mesh \"{0}\", "
758 "but no access region was defined although this is necessary for parallel runs. "
759 "Please define an access region using \"setMeshAccessRegion()\" before calling \"getMeshVertexSize()\".",
760 meshName);
761
762 auto result = mesh::countVerticesInBoundingBox(context.mesh, *context.userDefinedAccessRegion);
763
764 PRECICE_DEBUG("Filtered {} of {} vertices out on mesh {} due to the local access region. Mesh size in the access region: {}", context.mesh->nVertices() - result, context.mesh->nVertices(), meshName, result);
765 return result;
766 } else {
767 // For provided meshes and in case the api-access was not configured, we return here all vertices
768 PRECICE_WARN_IF(_accessor->isMeshReceived(meshName) && !_accessor->isDirectAccessAllowed(meshName),
769 "You are calling \"getMeshVertexSize()\" on a received mesh without api-access enabled (<receive-mesh name=\"{0}\" ... api-access=\"false\"/>). "
770 "Note that enabling api-access is required for this function to work properly with direct mesh access and just-in-time mappings.",
771 meshName);
772 return context.mesh->nVertices();
773 }
774}
775
778 std::string_view meshName)
779{
781 PRECICE_CHECK(_allowsRemeshing, "Cannot reset meshes. This feature needs to be enabled using <precice-configuration experimental=\"1\" allow-remeshing=\"1\">.");
782 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before resetMesh().");
783 PRECICE_TRACE(meshName);
785 PRECICE_CHECK(_couplingScheme->isCouplingOngoing(), "Cannot remesh after the last time window has been completed.");
786 PRECICE_CHECK(_couplingScheme->isTimeWindowComplete(), "Cannot remesh while subcycling or iterating. Remeshing is only allowed when the time window is completed.");
787 impl::MeshContext &context = _accessor->usedMeshContext(meshName);
788
789 PRECICE_DEBUG("Clear mesh positions for mesh \"{}\"", context.mesh->getName());
790 _meshLock.unlock(meshName);
791 context.mesh->clear();
792}
793
795 std::string_view meshName,
797{
798 PRECICE_TRACE(meshName);
800 MeshContext &context = _accessor->usedMeshContext(meshName);
801 auto &mesh = *context.mesh;
802 PRECICE_CHECK(position.size() == static_cast<unsigned long>(mesh.getDimensions()),
803 "Cannot set vertex for mesh \"{}\". Expected {} position components but found {}.", meshName, mesh.getDimensions(), position.size());
804 Event e{fmt::format("setMeshVertex.{}", meshName), profiling::Fundamental};
805 auto index = mesh.createVertex(Eigen::Map<const Eigen::VectorXd>{position.data(), mesh.getDimensions()}).getID();
807
808 const auto newSize = mesh.nVertices();
809 for (auto &context : _accessor->writeDataContexts()) {
810 if (context.getMeshName() == mesh.getName()) {
811 context.resizeBufferTo(newSize);
812 }
813 }
814
815 return index;
816}
817
819 std::string_view meshName,
822{
823 PRECICE_TRACE(meshName, positions.size(), ids.size());
825 MeshContext &context = _accessor->usedMeshContext(meshName);
826 auto &mesh = *context.mesh;
827
828 const auto meshDims = mesh.getDimensions();
829 const auto expectedPositionSize = ids.size() * meshDims;
830 PRECICE_CHECK(positions.size() == expectedPositionSize,
831 "Input sizes are inconsistent attempting to set vertices on {}D mesh \"{}\". "
832 "You passed {} vertex indices and {} position components, but we expected {} position components ({} x {}).",
833 meshDims, meshName, ids.size(), positions.size(), expectedPositionSize, ids.size(), meshDims);
834
835 Event e{fmt::format("setMeshVertices.{}", meshName), profiling::Fundamental};
836 const Eigen::Map<const Eigen::MatrixXd> posMatrix{
837 positions.data(), mesh.getDimensions(), static_cast<EIGEN_DEFAULT_DENSE_INDEX_TYPE>(ids.size())};
838 for (unsigned long i = 0; i < ids.size(); ++i) {
839 ids[i] = mesh.createVertex(posMatrix.col(i)).getID();
840 }
842
843 const auto newSize = mesh.nVertices();
844 for (auto &context : _accessor->writeDataContexts()) {
845 if (context.getMeshName() == mesh.getName()) {
846 context.resizeBufferTo(newSize);
847 }
848 }
849}
850
852 std::string_view meshName,
853 VertexID first,
854 VertexID second)
855{
856 PRECICE_TRACE(meshName, first, second);
858 MeshContext &context = _accessor->usedMeshContext(meshName);
860 mesh::PtrMesh &mesh = context.mesh;
862 PRECICE_CHECK(mesh->isValidVertexID(first), errorInvalidVertexID(first));
863 PRECICE_CHECK(mesh->isValidVertexID(second), errorInvalidVertexID(second));
864 Event e{fmt::format("setMeshEdge.{}", meshName), profiling::Fundamental};
865 mesh::Vertex &v0 = mesh->vertex(first);
866 mesh::Vertex &v1 = mesh->vertex(second);
867 mesh->createEdge(v0, v1);
868 }
869}
870
872 std::string_view meshName,
874{
875 PRECICE_TRACE(meshName, vertices.size());
877 MeshContext &context = _accessor->usedMeshContext(meshName);
879 return;
880 }
881
882 mesh::PtrMesh &mesh = context.mesh;
883 PRECICE_CHECK(vertices.size() % 2 == 0,
884 "Cannot interpret passed vertex IDs attempting to set edges of mesh \"{}\" . "
885 "You passed {} vertex indices, but we expected an even number.",
886 meshName, vertices.size());
887 {
888 auto end = vertices.end();
889 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
890 return !mesh->isValidVertexID(vid);
891 });
892 PRECICE_CHECK(first == end,
894 std::distance(vertices.begin(), first),
895 std::distance(vertices.begin(), last));
896 }
897
898 Event e{fmt::format("setMeshEdges.{}", meshName), profiling::Fundamental};
899
900 for (unsigned long i = 0; i < vertices.size() / 2; ++i) {
901 auto aid = vertices[2 * i];
902 auto bid = vertices[2 * i + 1];
903 mesh->createEdge(mesh->vertex(aid), mesh->vertex(bid));
904 }
905}
906
908 std::string_view meshName,
909 VertexID first,
910 VertexID second,
911 VertexID third)
912{
913 PRECICE_TRACE(meshName, first,
914 second, third);
915
917 MeshContext &context = _accessor->usedMeshContext(meshName);
919 mesh::PtrMesh &mesh = context.mesh;
921 PRECICE_CHECK(mesh->isValidVertexID(first), errorInvalidVertexID(first));
922 PRECICE_CHECK(mesh->isValidVertexID(second), errorInvalidVertexID(second));
923 PRECICE_CHECK(mesh->isValidVertexID(third), errorInvalidVertexID(third));
925 "setMeshTriangle() was called with repeated Vertex IDs ({}, {}, {}).",
926 first, second, third);
927 mesh::Vertex *vertices[3];
928 vertices[0] = &mesh->vertex(first);
929 vertices[1] = &mesh->vertex(second);
930 vertices[2] = &mesh->vertex(third);
932 vertices[1]->getCoords(), vertices[2]->getCoords())),
933 "setMeshTriangle() was called with vertices located at identical coordinates (IDs: {}, {}, {}).",
934 first, second, third);
935 Event e{fmt::format("setMeshTriangle.{}", meshName), profiling::Fundamental};
936 mesh::Edge *edges[3];
937 edges[0] = &mesh->createEdge(*vertices[0], *vertices[1]);
938 edges[1] = &mesh->createEdge(*vertices[1], *vertices[2]);
939 edges[2] = &mesh->createEdge(*vertices[2], *vertices[0]);
940
941 mesh->createTriangle(*edges[0], *edges[1], *edges[2]);
942 }
943}
944
946 std::string_view meshName,
948{
949 PRECICE_TRACE(meshName, vertices.size());
951 MeshContext &context = _accessor->usedMeshContext(meshName);
953 return;
954 }
955
956 mesh::PtrMesh &mesh = context.mesh;
957 PRECICE_CHECK(vertices.size() % 3 == 0,
958 "Cannot interpret passed vertex IDs attempting to set triangles of mesh \"{}\" . "
959 "You passed {} vertex indices, which isn't dividable by 3.",
960 meshName, vertices.size());
961 {
962 auto end = vertices.end();
963 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
964 return !mesh->isValidVertexID(vid);
965 });
966 PRECICE_CHECK(first == end,
968 std::distance(vertices.begin(), first),
969 std::distance(vertices.begin(), last));
970 }
971
972 Event e{fmt::format("setMeshTriangles.{}", meshName), profiling::Fundamental};
973
974 for (unsigned long i = 0; i < vertices.size() / 3; ++i) {
975 auto aid = vertices[3 * i];
976 auto bid = vertices[3 * i + 1];
977 auto cid = vertices[3 * i + 2];
978 mesh->createTriangle(mesh->vertex(aid),
979 mesh->vertex(bid),
980 mesh->vertex(cid));
981 }
982}
983
985 std::string_view meshName,
986 VertexID first,
987 VertexID second,
988 VertexID third,
989 VertexID fourth)
990{
991 PRECICE_TRACE(meshName, first,
992 second, third, fourth);
994 MeshContext &context = _accessor->usedMeshContext(meshName);
995 PRECICE_CHECK(context.mesh->getDimensions() == 3, "setMeshQuad is only possible for 3D meshes."
996 " Please set the mesh dimension to 3 in the preCICE configuration file.");
998 PRECICE_ASSERT(context.mesh);
999 mesh::Mesh &mesh = *(context.mesh);
1002 PRECICE_CHECK(mesh.isValidVertexID(second), errorInvalidVertexID(second));
1004 PRECICE_CHECK(mesh.isValidVertexID(fourth), errorInvalidVertexID(fourth));
1005
1006 auto vertexIDs = utils::make_array(first, second, third, fourth);
1007 PRECICE_CHECK(utils::unique_elements(vertexIDs), "The four vertex ID's are not unique. Please check that the vertices that form the quad are correct.");
1008
1009 auto coords = mesh::coordsFor(mesh, vertexIDs);
1011 "The four vertices that form the quad are not unique. The resulting shape may be a point, line or triangle."
1012 "Please check that the adapter sends the four unique vertices that form the quad, or that the mesh on the interface is composed of quads.");
1013
1014 auto convexity = math::geometry::isConvexQuad(coords);
1015 PRECICE_CHECK(convexity.convex, "The given quad is not convex. "
1016 "Please check that the adapter send the four correct vertices or that the interface is composed of quads.");
1017 auto reordered = utils::reorder_array(convexity.vertexOrder, mesh::vertexPtrsFor(mesh, vertexIDs));
1018
1019 Event e{fmt::format("setMeshQuad.{}", meshName), profiling::Fundamental};
1020
1021 // Vertices are now in the order: V0-V1-V2-V3-V0.
1022 // Use the shortest diagonal to split the quad into 2 triangles.
1023 // Vertices are now in V0-V1-V2-V3-V0 order. The new edge, e[4] is either 0-2 or 1-3
1024 double distance02 = (reordered[0]->getCoords() - reordered[2]->getCoords()).norm();
1025 double distance13 = (reordered[1]->getCoords() - reordered[3]->getCoords()).norm();
1026
1027 // The new edge, e[4], is the shortest diagonal of the quad
1028 if (distance02 <= distance13) {
1029 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[1]);
1030 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[3]);
1031 } else {
1032 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[0]);
1033 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[2]);
1034 }
1035 }
1036}
1037
1039 std::string_view meshName,
1041{
1042 PRECICE_TRACE(meshName, vertices.size());
1044 MeshContext &context = _accessor->usedMeshContext(meshName);
1046 return;
1047 }
1048
1049 mesh::Mesh &mesh = *(context.mesh);
1050 PRECICE_CHECK(vertices.size() % 4 == 0,
1051 "Cannot interpret passed vertex IDs attempting to set quads of mesh \"{}\" . "
1052 "You passed {} vertex indices, which isn't dividable by 4.",
1053 meshName, vertices.size());
1054 {
1055 auto end = vertices.end();
1056 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
1057 return !mesh.isValidVertexID(vid);
1058 });
1059 PRECICE_CHECK(first == end,
1061 std::distance(vertices.begin(), first),
1062 std::distance(vertices.begin(), last));
1063 }
1064
1065 for (unsigned long i = 0; i < vertices.size() / 4; ++i) {
1066 auto aid = vertices[4 * i];
1067 auto bid = vertices[4 * i + 1];
1068 auto cid = vertices[4 * i + 2];
1069 auto did = vertices[4 * i + 3];
1070
1071 auto vertexIDs = utils::make_array(aid, bid, cid, did);
1072 PRECICE_CHECK(utils::unique_elements(vertexIDs), "The four vertex ID's of the quad nr {} are not unique. Please check that the vertices that form the quad are correct.", i);
1073
1074 auto coords = mesh::coordsFor(mesh, vertexIDs);
1076 "The four vertices that form the quad nr {} are not unique. The resulting shape may be a point, line or triangle."
1077 "Please check that the adapter sends the four unique vertices that form the quad, or that the mesh on the interface is composed of quads.",
1078 i);
1079
1080 auto convexity = math::geometry::isConvexQuad(coords);
1081 PRECICE_CHECK(convexity.convex, "The given quad nr {} is not convex. "
1082 "Please check that the adapter send the four correct vertices or that the interface is composed of quads.",
1083 i);
1084 auto reordered = utils::reorder_array(convexity.vertexOrder, mesh::vertexPtrsFor(mesh, vertexIDs));
1085
1086 Event e{fmt::format("setMeshQuads.{}", meshName), profiling::Fundamental};
1087
1088 // Use the shortest diagonal to split the quad into 2 triangles.
1089 // Vertices are now in V0-V1-V2-V3-V0 order. The new edge, e[4] is either 0-2 or 1-3
1090 double distance02 = (reordered[0]->getCoords() - reordered[2]->getCoords()).norm();
1091 double distance13 = (reordered[1]->getCoords() - reordered[3]->getCoords()).norm();
1092
1093 if (distance02 <= distance13) {
1094 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[1]);
1095 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[3]);
1096 } else {
1097 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[0]);
1098 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[2]);
1099 }
1100 }
1101}
1102
1104 std::string_view meshName,
1105 VertexID first,
1106 VertexID second,
1107 VertexID third,
1108 VertexID fourth)
1109{
1110 PRECICE_TRACE(meshName, first, second, third, fourth);
1112 MeshContext &context = _accessor->usedMeshContext(meshName);
1113 PRECICE_CHECK(context.mesh->getDimensions() == 3, "setMeshTetrahedron is only possible for 3D meshes."
1114 " Please set the mesh dimension to 3 in the preCICE configuration file.");
1115 Event e{fmt::format("setMeshTetrahedron.{}", meshName), profiling::Fundamental};
1116
1118 mesh::PtrMesh &mesh = context.mesh;
1120 PRECICE_CHECK(mesh->isValidVertexID(first), errorInvalidVertexID(first));
1121 PRECICE_CHECK(mesh->isValidVertexID(second), errorInvalidVertexID(second));
1122 PRECICE_CHECK(mesh->isValidVertexID(third), errorInvalidVertexID(third));
1123 PRECICE_CHECK(mesh->isValidVertexID(fourth), errorInvalidVertexID(fourth));
1124 mesh::Vertex &A = mesh->vertex(first);
1125 mesh::Vertex &B = mesh->vertex(second);
1126 mesh::Vertex &C = mesh->vertex(third);
1127 mesh::Vertex &D = mesh->vertex(fourth);
1128
1129 mesh->createTetrahedron(A, B, C, D);
1130 }
1131}
1132
1134 std::string_view meshName,
1136{
1137 PRECICE_TRACE(meshName, vertices.size());
1139 MeshContext &context = _accessor->usedMeshContext(meshName);
1141 return;
1142 }
1143
1144 mesh::PtrMesh &mesh = context.mesh;
1145 PRECICE_CHECK(vertices.size() % 4 == 0,
1146 "Cannot interpret passed vertex IDs attempting to set quads of mesh \"{}\" . "
1147 "You passed {} vertex indices, which isn't dividable by 4.",
1148 meshName, vertices.size());
1149 {
1150 auto end = vertices.end();
1151 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
1152 return !mesh->isValidVertexID(vid);
1153 });
1154 PRECICE_CHECK(first == end,
1156 std::distance(vertices.begin(), first),
1157 std::distance(vertices.begin(), last));
1158 }
1159
1160 Event e{fmt::format("setMeshTetrahedra.{}", meshName), profiling::Fundamental};
1161
1162 for (unsigned long i = 0; i < vertices.size() / 4; ++i) {
1163 auto aid = vertices[4 * i];
1164 auto bid = vertices[4 * i + 1];
1165 auto cid = vertices[4 * i + 2];
1166 auto did = vertices[4 * i + 3];
1167 mesh->createTetrahedron(mesh->vertex(aid),
1168 mesh->vertex(bid),
1169 mesh->vertex(cid),
1170 mesh->vertex(did));
1171 }
1172}
1173
1175 std::string_view meshName,
1176 std::string_view dataName,
1179{
1180 PRECICE_TRACE(meshName, dataName, vertices.size());
1181 PRECICE_CHECK(_state != State::Finalized, "writeData(...) cannot be called after finalize().");
1182 PRECICE_CHECK(_state == State::Constructed || (_state == State::Initialized && isCouplingOngoing()), "Calling writeData(...) is forbidden if coupling is not ongoing, because the data you are trying to write will not be used anymore. You can fix this by always calling writeData(...) before the advance(...) call in your simulation loop or by using Participant::isCouplingOngoing() to implement a safeguard.");
1183 PRECICE_REQUIRE_DATA_WRITE(meshName, dataName);
1184 // Inconsistent sizes will be handled below
1185 if (vertices.empty() && values.empty()) {
1186 return;
1187 }
1188
1189 WriteDataContext &context = _accessor->writeDataContext(meshName, dataName);
1190
1191 const auto dataDims = context.getDataDimensions();
1192 const auto expectedDataSize = vertices.size() * dataDims;
1193 PRECICE_CHECK(expectedDataSize == values.size(),
1194 "Input sizes are inconsistent attempting to write {}D data \"{}\" to mesh \"{}\". "
1195 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1196 dataDims, dataName, meshName,
1197 vertices.size(), values.size(), expectedDataSize, dataDims, vertices.size());
1198
1199 // Sizes are correct at this point
1200 PRECICE_VALIDATE_DATA(values.data(), values.size()); // TODO Only take span
1201
1202 if (auto index = context.locateInvalidVertexID(vertices); index) {
1203 PRECICE_ERROR("Cannot write data \"{}\" to mesh \"{}\" due to invalid Vertex ID at vertices[{}]. "
1204 "Please make sure you only use the results from calls to setMeshVertex/Vertices().",
1205 dataName, meshName, *index);
1206 }
1207 Event e{fmt::format("writeData.{}_{}", meshName, dataName), profiling::Fundamental};
1208 context.writeValuesIntoDataBuffer(vertices, values);
1209}
1210
1212 std::string_view meshName,
1213 std::string_view dataName,
1215 double relativeReadTime,
1216 ::precice::span<double> values) const
1217{
1218 PRECICE_TRACE(meshName, dataName, vertices.size(), relativeReadTime);
1219 PRECICE_CHECK(_state != State::Constructed, "readData(...) cannot be called before initialize().");
1220 PRECICE_CHECK(_state != State::Finalized, "readData(...) cannot be called after finalize().");
1221 PRECICE_CHECK(math::smallerEquals(relativeReadTime, _couplingScheme->getNextTimeStepMaxSize()), "readData(...) cannot sample data outside of current time window.");
1222 PRECICE_CHECK(relativeReadTime >= 0, "readData(...) cannot sample data before the current time.");
1223 PRECICE_CHECK(isCouplingOngoing() || math::equals(relativeReadTime, 0.0), "Calling readData(...) with relativeReadTime = {} is forbidden if coupling is not ongoing. If coupling finished, only data for relativeReadTime = 0 is available. Please always use precice.getMaxTimeStepSize() to obtain the maximum allowed relativeReadTime.", relativeReadTime);
1224
1225 PRECICE_REQUIRE_DATA_READ(meshName, dataName);
1226
1227 PRECICE_CHECK(_meshLock.check(meshName),
1228 "Cannot read from mesh \"{}\" after it has been reset. Please read data before calling resetMesh().",
1229 meshName);
1230
1231 // Inconsistent sizes will be handled below
1232 if (vertices.empty() && values.empty()) {
1233 return;
1234 }
1235
1236 ReadDataContext &context = _accessor->readDataContext(meshName, dataName);
1237 PRECICE_CHECK(context.hasSamples(), "Data \"{}\" cannot be read from mesh \"{}\" as it contains no samples. "
1238 "This is typically a configuration issue of the data flow. "
1239 "Check if the data is correctly exchanged to this participant \"{}\" and mapped to mesh \"{}\".",
1240 dataName, meshName, _accessorName, meshName);
1241
1242 const auto dataDims = context.getDataDimensions();
1243 const auto expectedDataSize = vertices.size() * dataDims;
1244 PRECICE_CHECK(expectedDataSize == values.size(),
1245 "Input/Output sizes are inconsistent attempting to read {}D data \"{}\" from mesh \"{}\". "
1246 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1247 dataDims, dataName, meshName,
1248 vertices.size(), values.size(), expectedDataSize, dataDims, vertices.size());
1249
1250 if (auto index = context.locateInvalidVertexID(vertices); index) {
1251 PRECICE_ERROR("Cannot read data \"{}\" from mesh \"{}\" due to invalid Vertex ID at vertices[{}]. "
1252 "Please make sure you only use the results from calls to setMeshVertex/Vertices().",
1253 dataName, meshName, *index);
1254 }
1255
1256 Event e{fmt::format("readData.{}_{}", meshName, dataName), profiling::Fundamental};
1257
1258 double readTime = _couplingScheme->getTime() + relativeReadTime;
1259 context.readValues(vertices, readTime, values);
1260}
1261
1263 std::string_view meshName,
1264 std::string_view dataName,
1266 double relativeReadTime,
1267 ::precice::span<double> values) const
1268{
1270 PRECICE_TRACE(meshName, dataName, coordinates.size(), relativeReadTime);
1271 PRECICE_CHECK(_state != State::Constructed, "mapAndReadData(...) cannot be called before initialize().");
1272 PRECICE_CHECK(_state != State::Finalized, "mapAndReadData(...) cannot be called after finalize().");
1273 PRECICE_CHECK(math::smallerEquals(relativeReadTime, _couplingScheme->getNextTimeStepMaxSize()), "readData(...) cannot sample data outside of current time window.");
1274 PRECICE_CHECK(relativeReadTime >= 0, "mapAndReadData(...) cannot sample data before the current time.");
1275 PRECICE_CHECK(isCouplingOngoing() || math::equals(relativeReadTime, 0.0), "Calling mapAndReadData(...) with relativeReadTime = {} is forbidden if coupling is not ongoing. If coupling finished, only data for relativeReadTime = 0 is available. Please always use precice.getMaxTimeStepSize() to obtain the maximum allowed relativeReadTime.", relativeReadTime);
1276
1277 PRECICE_REQUIRE_DATA_READ(meshName, dataName);
1278 PRECICE_VALIDATE_DATA(coordinates.begin(), coordinates.size());
1279
1280 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1281 "This participant attempteded to map and read data (via \"mapAndReadData\") from mesh \"{0}\", "
1282 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1283 "mapAndReadData({0}, ...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1284 meshName);
1285 // If an access region is required, we have to check its existence
1286 bool requiresBB = requiresUserDefinedAccessRegion(meshName);
1287 PRECICE_CHECK(!requiresBB || (requiresBB && _accessor->meshContext(meshName).userDefinedAccessRegion),
1288 "The function \"mapAndReadData\" was called on mesh \"{0}\", "
1289 "but no access region was defined although this is necessary for parallel runs. "
1290 "Please define an access region using \"setMeshAccessRegion()\" before calling \"mapAndReadData()\".",
1291 meshName);
1292
1293 PRECICE_CHECK(!_accessor->meshContext(meshName).mesh->empty(), "This participant tries to mapAndRead data values for data \"{0}\" on mesh \"{1}\", but the mesh \"{1}\" is empty within the defined access region on this rank. "
1294 "How should the provided data values be read? Please make sure the mesh \"{1}\" is non-empty within the access region.",
1295 dataName, meshName);
1296
1297 // Inconsistent sizes will be handled below
1298 if (coordinates.empty() && values.empty()) {
1299 return;
1300 }
1301
1302 Event e{fmt::format("mapAndReadData.{}_{}", meshName, dataName), profiling::Fundamental};
1303
1304 // Note that meshName refers to a remote mesh
1305 ReadDataContext &dataContext = _accessor->readDataContext(meshName, dataName);
1306 const auto dataDims = dataContext.getDataDimensions();
1307 const auto dim = dataContext.getSpatialDimensions();
1308 const auto nVertices = (coordinates.size() / dim);
1309 MeshContext &context = _accessor->meshContext(meshName);
1310
1311 // Check that the vertex is actually within the defined access region
1312 context.checkVerticesInsideAccessRegion(coordinates, dim, "mapAndReadData");
1313
1314 // Make use of the read data context
1315 PRECICE_CHECK(nVertices * dataDims == values.size(),
1316 "Input sizes are inconsistent attempting to mapAndRead {}D data \"{}\" from mesh \"{}\". "
1317 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1318 dataDims, dataName, meshName,
1319 nVertices, values.size(), nVertices * dataDims, dataDims, nVertices);
1320
1321 double readTime = _couplingScheme->getTime() + relativeReadTime;
1322 dataContext.mapAndReadValues(coordinates, readTime, values);
1323}
1324
1326 std::string_view meshName,
1327 std::string_view dataName,
1330{
1332 PRECICE_TRACE(meshName, dataName, coordinates.size());
1333 PRECICE_CHECK(_state != State::Finalized, "writeAndMapData(...) cannot be called after finalize().");
1334 PRECICE_CHECK(_state == State::Constructed || (_state == State::Initialized && isCouplingOngoing()), "Calling writeAndMapData(...) is forbidden if coupling is not ongoing, because the data you are trying to write will not be used anymore. You can fix this by always calling writeAndMapData(...) before the advance(...) call in your simulation loop or by using Participant::isCouplingOngoing() to implement a safeguard.");
1335 PRECICE_REQUIRE_DATA_WRITE(meshName, dataName);
1336
1337 PRECICE_VALIDATE_DATA(coordinates.begin(), coordinates.size());
1338 PRECICE_VALIDATE_DATA(values.data(), values.size());
1339 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1340 "This participant attempteded to map and read data (via \"writeAndMapData\") from mesh \"{0}\", "
1341 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1342 "writeAndMapData({0}, ...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1343 meshName);
1344 // If an access region is required, we have to check its existence
1345 bool requiresBB = requiresUserDefinedAccessRegion(meshName);
1346 PRECICE_CHECK(!requiresBB || (requiresBB && _accessor->meshContext(meshName).userDefinedAccessRegion),
1347 "The function \"writeAndMapData\" was called on mesh \"{0}\", "
1348 "but no access region was defined although this is necessary for parallel runs. "
1349 "Please define an access region using \"setMeshAccessRegion()\" before calling \"writeAndMapData()\".",
1350 meshName);
1351
1352 // Inconsistent sizes will be handled below
1353 if (coordinates.empty() && values.empty()) {
1354 return;
1355 }
1356
1357 Event e{fmt::format("writeAndMapData.{}_{}", meshName, dataName), profiling::Fundamental};
1358
1359 // Note that meshName refers here typically to a remote mesh
1360 WriteDataContext &dataContext = _accessor->writeDataContext(meshName, dataName);
1361 const auto dataDims = dataContext.getDataDimensions();
1362 const auto dim = dataContext.getSpatialDimensions();
1363 const auto nVertices = (coordinates.size() / dim);
1364 MeshContext &context = _accessor->meshContext(meshName);
1365
1366 // Check that the vertex is actually within the defined access region
1367 context.checkVerticesInsideAccessRegion(coordinates, dim, "writeAndMapData");
1368
1369 PRECICE_CHECK(nVertices * dataDims == values.size(),
1370 "Input sizes are inconsistent attempting to write {}D data \"{}\" to mesh \"{}\". "
1371 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1372 dataDims, dataName, meshName,
1373 nVertices, values.size(), nVertices * dataDims, dataDims, nVertices);
1374
1375 PRECICE_CHECK(!context.mesh->empty(), "This participant tries to mapAndWrite data values for data \"{0}\" on mesh \"{1}\", but the mesh \"{1}\" is empty within the defined access region on this rank. "
1376 "Where should the provided data go? Please make sure the mesh \"{1}\" is non-empty within the access region.",
1377 dataName, meshName);
1378 dataContext.writeAndMapValues(coordinates, values);
1379}
1380
1382 std::string_view meshName,
1383 std::string_view dataName,
1386{
1388
1389 // Asserts and checks
1390 PRECICE_TRACE(meshName, dataName, vertices.size());
1391 PRECICE_CHECK(_state != State::Finalized, "writeGradientData(...) cannot be called after finalize().");
1392 PRECICE_REQUIRE_DATA_WRITE(meshName, dataName);
1393
1394 // Inconsistent sizes will be handled below
1395 if ((vertices.empty() && gradients.empty()) || !requiresGradientDataFor(meshName, dataName)) {
1396 return;
1397 }
1398
1399 // Get the data
1400 WriteDataContext &context = _accessor->writeDataContext(meshName, dataName);
1401
1402 // Check if the Data object of given mesh has been initialized with gradient data
1403 PRECICE_CHECK(context.hasGradient(), "Data \"{}\" has no gradient values available. Please set the gradient flag to true under the data attribute in the configuration file.", dataName);
1404
1405 if (auto index = context.locateInvalidVertexID(vertices); index) {
1406 PRECICE_ERROR("Cannot write gradient data \"{}\" to mesh \"{}\" due to invalid Vertex ID at vertices[{}]. "
1407 "Please make sure you only use the results from calls to setMeshVertex/Vertices().",
1408 dataName, meshName, *index);
1409 }
1410
1411 const auto dataDims = context.getDataDimensions();
1412 const auto meshDims = context.getSpatialDimensions();
1413 const auto gradientComponents = meshDims * dataDims;
1414 const auto expectedComponents = vertices.size() * gradientComponents;
1415 PRECICE_CHECK(expectedComponents == gradients.size(),
1416 "Input sizes are inconsistent attempting to write gradient for data \"{}\" to mesh \"{}\". "
1417 "A single gradient/Jacobian for {}D data on a {}D mesh has {} components. "
1418 "You passed {} vertex indices and {} gradient components, but we expected {} gradient components. ",
1419 dataName, meshName,
1420 dataDims, meshDims, gradientComponents,
1421 vertices.size(), gradients.size(), expectedComponents);
1422
1423 PRECICE_VALIDATE_DATA(gradients.data(), gradients.size());
1424
1425 Event e{fmt::format("writeGradientData.{}_{}", meshName, dataName), profiling::Fundamental};
1426
1427 context.writeGradientsIntoDataBuffer(vertices, gradients);
1428}
1429
1431 const std::string_view meshName,
1432 ::precice::span<const double> boundingBox) const
1433{
1434 PRECICE_TRACE(meshName, boundingBox.size());
1435 PRECICE_REQUIRE_MESH_USE(meshName);
1436 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1437 "This participant attempteded to set an access region (via \"setMeshAccessRegion\") on mesh \"{0}\", "
1438 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1439 "setMeshAccessRegion(...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1440 meshName);
1441 PRECICE_CHECK(_state != State::Finalized, "setMeshAccessRegion() cannot be called after finalize().");
1442 PRECICE_CHECK(_state != State::Initialized, "setMeshAccessRegion() needs to be called before initialize().");
1443
1444 // Get the related mesh
1445 MeshContext &context = _accessor->meshContext(meshName);
1446
1447 PRECICE_CHECK(!context.userDefinedAccessRegion, "A mesh access region was already defined for mesh \"{}\". setMeshAccessRegion may only be called once per mesh.", context.mesh->getName());
1448 mesh::PtrMesh mesh(context.mesh);
1449 int dim = mesh->getDimensions();
1450 PRECICE_CHECK(boundingBox.size() == static_cast<unsigned long>(dim) * 2,
1451 "Incorrect amount of bounding box components attempting to set the bounding box of {}D mesh \"{}\" . "
1452 "You passed {} limits, but we expected {} ({}x2).",
1453 dim, meshName, boundingBox.size(), dim * 2, dim);
1454
1455 // Transform bounds into a suitable format
1456 PRECICE_DEBUG("Define bounding box");
1457 std::vector<double> bounds(dim * 2);
1458
1459 for (int d = 0; d < dim; ++d) {
1460 // Check that min is lower or equal to max
1461 PRECICE_CHECK(boundingBox[2 * d] <= boundingBox[2 * d + 1], "Your bounding box is ill defined, i.e. it has a negative volume. The required format is [x_min, x_max...]");
1462 bounds[2 * d] = boundingBox[2 * d];
1463 bounds[2 * d + 1] = boundingBox[2 * d + 1];
1464 }
1465 // Create a bounding box
1467 // Expand the mesh associated bounding box
1468 mesh->expandBoundingBox(*context.userDefinedAccessRegion.get());
1469}
1470
1472 const std::string_view meshName,
1474 ::precice::span<double> coordinates) const
1475{
1476 PRECICE_TRACE(meshName, ids.size(), coordinates.size());
1477 PRECICE_REQUIRE_MESH_USE(meshName);
1478 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1479 "This participant attempteded to get mesh vertex IDs and coordinates (via \"getMeshVertexIDsAndCoordinates\") from mesh \"{0}\", "
1480 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1481 "getMeshVertexIDsAndCoordinates(...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1482 meshName);
1483 // If an access region is required, we have to check its existence
1484 bool requiresBB = requiresUserDefinedAccessRegion(meshName);
1485 PRECICE_CHECK(!requiresBB || (requiresBB && _accessor->meshContext(meshName).userDefinedAccessRegion),
1486 "The function \"getMeshVertexIDsAndCoordinates\" was called on mesh \"{0}\", "
1487 "but no access region was defined although this is necessary for parallel runs. "
1488 "Please define an access region using \"setMeshAccessRegion()\" before calling \"getMeshVertexIDsAndCoordinates()\".",
1489 meshName);
1490
1491 PRECICE_DEBUG("Get {} mesh vertices with IDs", ids.size());
1492
1493 // Check, if the requested mesh data has already been received. Otherwise, the function call doesn't make any sense
1494 PRECICE_CHECK((_state == State::Initialized) || _accessor->isMeshProvided(meshName), "initialize() has to be called before accessing"
1495 " data of the received mesh \"{}\" on participant \"{}\".",
1496 meshName, _accessor->getName());
1497
1498 if (ids.empty() && coordinates.empty()) {
1499 return;
1500 }
1501
1502 Event e{fmt::format("getMeshVertexIDsAndCoordinates.{}", meshName), profiling::Fundamental};
1503
1504 const MeshContext &context = _accessor->meshContext(meshName);
1505
1506 auto filteredVertices = context.filterVerticesToLocalAccessRegion(requiresBB);
1507 const auto meshSize = filteredVertices.size();
1508
1509 const mesh::PtrMesh mesh(context.mesh);
1510 const auto meshDims = mesh->getDimensions();
1511 PRECICE_CHECK(ids.size() == meshSize,
1512 "Output size is incorrect attempting to get vertex ids of {}D mesh \"{}\". "
1513 "You passed {} vertex indices, but we expected {}. "
1514 "Use getMeshVertexSize(\"{}\") to receive the required amount of vertices.",
1515 meshDims, meshName, ids.size(), meshSize, meshName);
1516 const auto expectedCoordinatesSize = static_cast<unsigned long>(meshDims * meshSize);
1517 PRECICE_CHECK(coordinates.size() == expectedCoordinatesSize,
1518 "Output size is incorrect attempting to get vertex coordinates of {}D mesh \"{}\". "
1519 "You passed {} coordinate components, but we expected {} ({}x{}). "
1520 "Use getMeshVertexSize(\"{}\") and getMeshDimensions(\"{}\") to receive the required amount components",
1521 meshDims, meshName, coordinates.size(), expectedCoordinatesSize, meshSize, meshDims, meshName, meshName);
1522
1523 PRECICE_ASSERT(ids.size() <= mesh->nVertices(), "The queried size exceeds the number of available points.");
1524
1525 Eigen::Map<Eigen::MatrixXd> posMatrix{
1526 coordinates.data(), mesh->getDimensions(), static_cast<EIGEN_DEFAULT_DENSE_INDEX_TYPE>(ids.size())};
1527
1528 for (unsigned long i = 0; i < ids.size(); i++) {
1529 auto localID = filteredVertices[i].get().getID();
1530 PRECICE_ASSERT(mesh->isValidVertexID(localID), i, localID);
1531 ids[i] = localID;
1532 posMatrix.col(i) = filteredVertices[i].get().getCoords();
1533 }
1534}
1535
1538{
1539 PRECICE_TRACE();
1540 for (const auto &m2nConf : config->m2ns()) {
1541 if (m2nConf.acceptor != _accessorName && m2nConf.connector != _accessorName) {
1542 continue;
1543 }
1544
1545 std::string comPartner("");
1546 bool isRequesting;
1547 if (m2nConf.acceptor == _accessorName) {
1548 comPartner = m2nConf.connector;
1549 isRequesting = true;
1550 } else {
1551 comPartner = m2nConf.acceptor;
1552 isRequesting = false;
1553 }
1554
1555 PRECICE_ASSERT(!comPartner.empty());
1556 for (const impl::PtrParticipant &participant : _participants) {
1557 if (participant->getName() == comPartner) {
1558 PRECICE_ASSERT(not utils::contained(comPartner, _m2ns), comPartner);
1559 PRECICE_ASSERT(m2nConf.m2n);
1560
1561 _m2ns[comPartner] = [&] {
1562 m2n::BoundM2N bound;
1563 bound.m2n = m2nConf.m2n;
1564 bound.localName = _accessorName;
1565 bound.remoteName = comPartner;
1566 bound.isRequesting = isRequesting;
1567 return bound;
1568 }();
1569 }
1570 }
1571 }
1572}
1573
1575 const m2n::M2NConfiguration::SharedPointer &m2nConfig)
1576{
1577 PRECICE_TRACE();
1578 for (MeshContext *context : _accessor->usedMeshContexts()) {
1579
1580 if (context->provideMesh) { // Accessor provides mesh
1581 PRECICE_CHECK(context->receiveMeshFrom.empty(),
1582 "Participant \"{}\" cannot provide and receive mesh {}!",
1583 _accessorName, context->mesh->getName());
1584
1585 context->partition = partition::PtrPartition(new partition::ProvidedPartition(context->mesh));
1586
1587 for (auto &receiver : _participants) {
1588 for (auto &receiverContext : receiver->usedMeshContexts()) {
1589 if (receiverContext->receiveMeshFrom == _accessorName && receiverContext->mesh->getName() == context->mesh->getName()) {
1590 // meshRequirement has to be copied from "from" to provide", since
1591 // mapping are only defined at "provide"
1592 if (receiverContext->meshRequirement > context->meshRequirement) {
1593 context->meshRequirement = receiverContext->meshRequirement;
1594 }
1595
1596 m2n::PtrM2N m2n = m2nConfig->getM2N(receiver->getName(), _accessorName);
1597 m2n->createDistributedCommunication(context->mesh);
1598 context->partition->addM2N(m2n);
1599 }
1600 }
1601 }
1602
1603 } else { // Accessor receives mesh
1604 std::string receiver(_accessorName);
1605 std::string provider(context->receiveMeshFrom);
1606
1607 PRECICE_DEBUG("Receiving mesh from {}", provider);
1608
1609 context->partition = partition::PtrPartition(new partition::ReceivedPartition(context->mesh, context->geoFilter, context->safetyFactor, context->allowDirectAccess));
1610
1611 m2n::PtrM2N m2n = m2nConfig->getM2N(receiver, provider);
1612 m2n->createDistributedCommunication(context->mesh);
1613 context->partition->addM2N(m2n);
1614 for (const MappingContext &mappingContext : context->fromMappingContexts) {
1615 context->partition->addFromMapping(mappingContext.mapping);
1616 }
1617 for (const MappingContext &mappingContext : context->toMappingContexts) {
1618 context->partition->addToMapping(mappingContext.mapping);
1619 }
1620 }
1621 }
1622}
1623
1625{
1626 // sort meshContexts by name, for communication in right order.
1627 std::sort(_accessor->usedMeshContexts().begin(), _accessor->usedMeshContexts().end(),
1628 [](MeshContext const *const lhs, MeshContext const *const rhs) -> bool {
1629 return lhs->mesh->getName() < rhs->mesh->getName();
1630 });
1631
1632 for (MeshContext *meshContext : _accessor->usedMeshContexts()) {
1633 if (meshContext->provideMesh) // provided meshes need their bounding boxes already for the re-partitioning
1634 meshContext->mesh->computeBoundingBox();
1635
1636 meshContext->clearMappings();
1637 }
1638
1639 for (MeshContext *meshContext : _accessor->usedMeshContexts()) {
1640 meshContext->partition->compareBoundingBoxes();
1641 }
1642}
1643
1645{
1646 // We need to do this in two loops: First, communicate the mesh and later compute the partition.
1647 // Originally, this was done in one loop. This however gave deadlock if two meshes needed to be communicated cross-wise.
1648 // Both loops need a different sorting
1649
1650 auto &contexts = _accessor->usedMeshContexts();
1651
1652 std::sort(contexts.begin(), contexts.end(),
1653 [](MeshContext const *const lhs, MeshContext const *const rhs) -> bool {
1654 return lhs->mesh->getName() < rhs->mesh->getName();
1655 });
1656
1657 for (MeshContext *meshContext : contexts) {
1658 meshContext->partition->communicate();
1659 }
1660
1661 // for two-level initialization, there is also still communication in partition::compute()
1662 // therefore, we cannot resort here.
1663 // @todo this hacky solution should be removed as part of #633
1664 bool resort = true;
1665 for (auto &m2nPair : _m2ns) {
1666 if (m2nPair.second.m2n->usesTwoLevelInitialization()) {
1667 resort = false;
1668 break;
1669 }
1670 }
1671
1672 if (resort) {
1673 // pull provided meshes up front, to have them ready for the decomposition of the received meshes (for the mappings)
1674 std::stable_partition(contexts.begin(), contexts.end(),
1675 [](MeshContext const *const meshContext) -> bool {
1676 return meshContext->provideMesh;
1677 });
1678 }
1679
1680 for (MeshContext *meshContext : contexts) {
1681 meshContext->partition->compute();
1682 if (not meshContext->provideMesh) { // received mesh can only compute their bounding boxes here
1683 meshContext->mesh->computeBoundingBox();
1684 }
1685
1686 meshContext->mesh->allocateDataValues();
1687
1688 const auto requiredSize = meshContext->mesh->nVertices();
1689 for (auto &context : _accessor->writeDataContexts()) {
1690 if (context.getMeshName() == meshContext->mesh->getName()) {
1691 context.resizeBufferTo(requiredSize);
1692 }
1693 }
1694 }
1695}
1696
1698{
1699 PRECICE_TRACE();
1700 bool anyMappingChanged = false;
1701 for (impl::MappingContext &context : contexts) {
1702 if (not context.mapping->hasComputedMapping()) {
1703 PRECICE_INFO_IF(context.configuredWithAliasTag,
1704 "Automatic RBF mapping alias from mesh \"{}\" to mesh \"{}\" in \"{}\" direction resolves to \"{}\" .",
1705 context.mapping->getInputMesh()->getName(), context.mapping->getOutputMesh()->getName(), mappingType, context.mapping->getName());
1706 PRECICE_INFO("Computing \"{}\" mapping from mesh \"{}\" to mesh \"{}\" in \"{}\" direction.",
1707 context.mapping->getName(), context.mapping->getInputMesh()->getName(), context.mapping->getOutputMesh()->getName(), mappingType);
1708 context.mapping->computeMapping();
1709 anyMappingChanged = true;
1710 }
1711 }
1712 if (anyMappingChanged) {
1713 _accessor->initializeMappingDataCache(mappingType);
1714 }
1715}
1716
1718{
1719 PRECICE_TRACE();
1720 if (!_accessor->hasWriteMappings()) {
1721 return;
1722 }
1723
1724 Event e("mapping", profiling::Fundamental);
1725 computeMappings(_accessor->writeMappingContexts(), "write");
1726 for (auto &context : _accessor->writeDataContexts()) {
1727 if (context.hasMapping()) {
1728 PRECICE_DEBUG("Map initial write data \"{}\" from mesh \"{}\"", context.getDataName(), context.getMeshName());
1729 _executedWriteMappings += context.mapData(std::nullopt, true);
1730 }
1731 }
1732}
1733
1735{
1736 PRECICE_TRACE();
1737 if (!_accessor->hasWriteMappings()) {
1738 return;
1739 }
1740
1741 Event e("mapping", profiling::Fundamental);
1742 computeMappings(_accessor->writeMappingContexts(), "write");
1743 for (auto &context : _accessor->writeDataContexts()) {
1744 if (context.hasMapping()) {
1745 PRECICE_DEBUG("Map write data \"{}\" from mesh \"{}\"", context.getDataName(), context.getMeshName());
1746 _executedWriteMappings += context.mapData(after);
1747 }
1748 }
1749}
1750
1751void ParticipantImpl::trimReadMappedData(double startOfTimeWindow, bool isTimeWindowComplete, const cplscheme::ImplicitData &fromData)
1752{
1753 PRECICE_TRACE();
1754 for (auto &context : _accessor->readDataContexts()) {
1755 if (context.hasMapping()) {
1757 // For serial implicit second, we need to discard everything before startOfTimeWindow to preserve the time window start
1758 // For serial implicit first, we need to discard everything as everything is new
1759 // For parallel implicit, we need to discard everything as everything is new
1760 context.clearToDataFor(fromData);
1761 } else {
1762 context.trimToDataAfterFor(fromData, startOfTimeWindow);
1763 }
1764 }
1765 }
1766}
1767
1769{
1770 PRECICE_TRACE();
1771 if (!_accessor->hasReadMappings()) {
1772 return;
1773 }
1774
1775 Event e("mapping", profiling::Fundamental);
1776 computeMappings(_accessor->readMappingContexts(), "read");
1777 for (auto &context : _accessor->readDataContexts()) {
1778 if (context.hasMapping()) {
1779 PRECICE_DEBUG("Map initial read data \"{}\" to mesh \"{}\"", context.getDataName(), context.getMeshName());
1780 // We always ensure that all read data was mapped
1781 _executedReadMappings += context.mapData(std::nullopt, true);
1782 }
1783 }
1784}
1785
1787{
1788 PRECICE_TRACE();
1789 if (!_accessor->hasReadMappings()) {
1790 return;
1791 }
1792
1793 Event e("mapping", profiling::Fundamental);
1794 computeMappings(_accessor->readMappingContexts(), "read");
1795 for (auto &context : _accessor->readDataContexts()) {
1796 if (context.hasMapping()) {
1797 PRECICE_DEBUG("Map read data \"{}\" to mesh \"{}\"", context.getDataName(), context.getMeshName());
1798 // We always ensure that all read data was mapped
1799 _executedReadMappings += context.mapData();
1800 }
1801 }
1802}
1803
1805{
1806 PRECICE_TRACE();
1807 for (action::PtrAction &action : _accessor->actions()) {
1808 if (timings.find(action->getTiming()) != timings.end()) {
1809 action->performAction();
1810 }
1811 }
1812}
1813
1815{
1816 PRECICE_TRACE();
1817 if (!_accessor->hasExports()) {
1818 return;
1819 }
1820 PRECICE_DEBUG("Handle exports");
1821 profiling::Event e{"handleExports"};
1822
1823 if (timing == ExportTiming::Initial) {
1824 _accessor->exportInitial();
1825 return;
1826 }
1827
1829 exp.timewindow = _couplingScheme->getTimeWindows() - 1;
1830 exp.iteration = _numberAdvanceCalls;
1831 exp.complete = _couplingScheme->isTimeWindowComplete();
1832 exp.final = !_couplingScheme->isCouplingOngoing();
1833 exp.time = _couplingScheme->getTime();
1834 _accessor->exportIntermediate(exp);
1835}
1836
1838{
1839 PRECICE_TRACE();
1840 for (auto &context : _accessor->writeDataContexts()) {
1841 // reset the buffered data here
1842 context.resetBufferedData();
1843 }
1844}
1845
1847 const config::Configuration &config)
1848{
1849 const auto &partConfig = config.getParticipantConfiguration();
1850 for (const PtrParticipant &participant : partConfig->getParticipants()) {
1851 if (participant->getName() == _accessorName) {
1852 return participant;
1853 }
1854 }
1855 PRECICE_ERROR("This participant's name, which was specified in the constructor of the preCICE interface as \"{}\", "
1856 "is not defined in the preCICE configuration. "
1857 "Please double-check the correct spelling.",
1859}
1860
1871
1872void ParticipantImpl::syncTimestep(double computedTimeStepSize)
1873{
1875 Event e("syncTimestep", profiling::Fundamental);
1877 utils::IntraComm::getCommunication()->send(computedTimeStepSize, 0);
1878 } else {
1880 for (Rank secondaryRank : utils::IntraComm::allSecondaryRanks()) {
1881 double dt;
1882 utils::IntraComm::getCommunication()->receive(dt, secondaryRank);
1883 PRECICE_CHECK(math::equals(dt, computedTimeStepSize),
1884 "Found ambiguous values for the time step size passed to preCICE in \"advance\". On rank {}, the value is {}, while on rank 0, the value is {}.",
1885 secondaryRank, dt, computedTimeStepSize);
1886 }
1887 }
1888}
1889
1891{
1892 PRECICE_DEBUG("Advance coupling scheme");
1893 Event e("advanceCoupling", profiling::Fundamental);
1894 // Orchestrate local and remote mesh changes
1895 std::vector<MeshID> localChanges;
1896
1897 [[maybe_unused]] auto remoteChanges1 = _couplingScheme->firstSynchronization(localChanges);
1898 _couplingScheme->firstExchange();
1899 // Orchestrate remote mesh changes (local ones were handled in the first sync)
1900 [[maybe_unused]] auto remoteChanges2 = _couplingScheme->secondSynchronization();
1901 _couplingScheme->secondExchange();
1902}
1903
1905{
1906 // Optionally apply some final ping-pong to sync solver that run e.g. with a uni-directional coupling
1907 // afterwards close connections
1908 PRECICE_INFO("{} {}communication channels",
1909 (_waitInFinalize ? "Synchronize participants and close" : "Close"),
1910 (close == CloseChannels::Distributed ? "distributed " : ""));
1911 std::string ping = "ping";
1912 std::string pong = "pong";
1913 for (auto &iter : _m2ns) {
1914 auto bm2n = iter.second;
1915 if (!bm2n.m2n->isConnected()) {
1916 PRECICE_DEBUG("Skipping closure of defective connection with {}", bm2n.remoteName);
1917 continue;
1918 }
1920 auto comm = bm2n.m2n->getPrimaryRankCommunication();
1921 PRECICE_DEBUG("Synchronizing primary rank with {}", bm2n.remoteName);
1922 if (bm2n.isRequesting) {
1923 comm->send(ping, 0);
1924 std::string receive = "init";
1925 comm->receive(receive, 0);
1926 PRECICE_ASSERT(receive == pong);
1927 } else {
1928 std::string receive = "init";
1929 comm->receive(receive, 0);
1930 PRECICE_ASSERT(receive == ping);
1931 comm->send(pong, 0);
1932 }
1933 }
1934 if (close == CloseChannels::Distributed) {
1935 PRECICE_DEBUG("Closing distributed communication with {}", bm2n.remoteName);
1936 bm2n.m2n->closeDistributedConnections();
1937 } else {
1938 PRECICE_DEBUG("Closing communication with {}", bm2n.remoteName);
1939 bm2n.m2n->closeConnection();
1940 }
1941 }
1942}
1943
1945{
1946 return _accessor->isMeshReceived(meshName) && utils::IntraComm::isParallel();
1947}
1948
1949const mesh::Mesh &ParticipantImpl::mesh(const std::string &meshName) const
1950{
1951 PRECICE_TRACE(meshName);
1952 return *_accessor->usedMeshContext(meshName).mesh;
1953}
1954
1962
1963// Reinitialization
1964
1966{
1967 PRECICE_TRACE();
1969 Event e("remesh.exchangeLocalMeshChanges", profiling::Synchronize);
1970
1971 // Gather local changes
1972 std::vector<double> localMeshChanges;
1973 for (auto context : _accessor->usedMeshContexts()) {
1974 localMeshChanges.push_back(_meshLock.check(context->mesh->getName()) ? 0.0 : 1.0);
1975 }
1976 PRECICE_DEBUG("Mesh changes of rank: {}", localMeshChanges);
1977
1978 // TODO implement int version of allreduceSum
1979 std::vector<double> totalMeshChanges(localMeshChanges.size(), 0.0);
1980 utils::IntraComm::allreduceSum(localMeshChanges, totalMeshChanges);
1981
1982 // Convert the doubles to int
1983 MeshChanges totalMeshChangesInt(totalMeshChanges.begin(), totalMeshChanges.end());
1984 PRECICE_DEBUG("Mesh changes of participant: {}", totalMeshChangesInt);
1985 return totalMeshChangesInt;
1986}
1987
1989{
1990 auto meshContexts = _accessor->usedMeshContexts();
1991 for (std::size_t i = 0; i < totalMeshChanges.size(); ++i) {
1992 if (totalMeshChanges[i] > 0.0) {
1993 meshContexts[i]->mesh->clearDataStamples();
1994 }
1995 }
1996}
1997
1998bool ParticipantImpl::reinitHandshake(bool requestReinit) const
1999{
2000 PRECICE_TRACE();
2002 Event e("remesh.exchangeRemoteMeshChanges", profiling::Synchronize);
2003
2005 PRECICE_DEBUG("Remeshing is{} required by this participant.", (requestReinit ? "" : " not"));
2006
2007 bool swarmReinitRequired = requestReinit;
2008 for (auto &iter : _m2ns) {
2009 PRECICE_DEBUG("Coordinating remeshing with {}", iter.first);
2010 bool received = false;
2011 auto &comm = *iter.second.m2n->getPrimaryRankCommunication();
2012 if (iter.second.isRequesting) {
2013 comm.send(requestReinit, 0);
2014 comm.receive(received, 0);
2015 } else {
2016 comm.receive(received, 0);
2017 comm.send(requestReinit, 0);
2018 }
2019 swarmReinitRequired |= received;
2020 }
2021 PRECICE_DEBUG("Coordinated that overall{} remeshing is required.", (swarmReinitRequired ? "" : " no"));
2022
2023 utils::IntraComm::broadcast(swarmReinitRequired);
2024 return swarmReinitRequired;
2025 } else {
2026 bool swarmReinitRequired = false;
2027 utils::IntraComm::broadcast(swarmReinitRequired);
2028 return swarmReinitRequired;
2029 }
2030}
2031
2033{
2034 PRECICE_CHECK(std::find(sectionName.begin(), sectionName.end(), '/') == sectionName.end(),
2035 "The provided section name \"{}\" may not contain a forward-slash \"/\"",
2036 sectionName);
2037 _userEvents.emplace_back(sectionName, profiling::Fundamental);
2038}
2039
2041{
2042 PRECICE_CHECK(!_userEvents.empty(), "There is no user-started event to stop.");
2043 _userEvents.pop_back();
2044}
2045
2046} // namespace precice::impl
#define PRECICE_ERROR(...)
Definition LogMacros.hpp:16
#define PRECICE_WARN_IF(condition,...)
Definition LogMacros.hpp:18
#define PRECICE_DEBUG(...)
Definition LogMacros.hpp:61
#define PRECICE_TRACE(...)
Definition LogMacros.hpp:92
#define PRECICE_INFO_IF(condition,...)
Definition LogMacros.hpp:25
#define PRECICE_INFO(...)
Definition LogMacros.hpp:14
#define PRECICE_CHECK(check,...)
Definition LogMacros.hpp:32
unsigned int index
std::string name
#define PRECICE_VALIDATE_DATA_NAME(mesh, data)
#define PRECICE_REQUIRE_DATA_WRITE(mesh, data)
#define PRECICE_REQUIRE_MESH_USE(name)
#define PRECICE_REQUIRE_DATA_READ(mesh, data)
#define PRECICE_REQUIRE_MESH_MODIFY(name)
#define PRECICE_VALIDATE_DATA(data, size)
#define PRECICE_EXPERIMENTAL_API()
#define PRECICE_VALIDATE_MESH_NAME(name)
T accumulate(T... args)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
T begin(T... args)
Main class for preCICE XML configuration tree.
const PtrParticipantConfiguration & getParticipantConfiguration() const
const mesh::PtrMeshConfiguration getMeshConfiguration() const
const m2n::M2NConfiguration::SharedPointer getM2NConfiguration() const
bool waitInFinalize() const
Returns whether participants wait for each other in finalize.
const cplscheme::PtrCouplingSchemeConfiguration getCouplingSchemeConfiguration() const
bool allowsExperimental() const
Returns whether experimental features are allowed or not.
bool allowsRemeshing() const
Returns whether experimental remeshing is allowed or not.
xml::XMLTag & getXMLTag()
Returns root xml tag to start the automatic configuration process.
@ WriteCheckpoint
Is the participant required to write a checkpoint?
@ ReadCheckpoint
Is the participant required to read a previously written checkpoint?
@ InitializeData
Is the initialization of coupling data required?
static void finalize()
Definition Ginkgo.cpp:36
bool hasGradient() const
Returns whether _providedData has gradient.
int getDataDimensions() const
Get the dimensions of _providedData.
int getSpatialDimensions() const
Get the spatial dimensions of _providedData.
std::optional< std::size_t > locateInvalidVertexID(const Container &c)
CloseChannels
Which channels to close in closeCommunicationChannels()
void writeGradientData(std::string_view meshName, std::string_view dataName, ::precice::span< const VertexID > vertices, ::precice::span< const double > gradients)
Writes vector gradient data to a mesh.
int getMeshVertexSize(std::string_view meshName) const
Returns the number of vertices of a mesh.
utils::MultiLock< std::string > _meshLock
void setMeshQuad(std::string_view meshName, VertexID first, VertexID second, VertexID third, VertexID fourth)
Sets a planar surface mesh quadrangle from vertex IDs.
void setMeshTetrahedra(std::string_view meshName, ::precice::span< const VertexID > vertices)
Sets multiple mesh tetrahedra from vertex IDs.
bool requiresGradientDataFor(std::string_view meshName, std::string_view dataName) const
Checks if the given data set requires gradient data. We check if the data object has been initialized...
int getDataDimensions(std::string_view meshName, std::string_view dataName) const
Returns the spatial dimensionality of the given data on the given mesh.
impl::PtrParticipant determineAccessingParticipant(const config::Configuration &config)
Determines participant accessing this interface from the configuration.
MeshChanges getTotalMeshChanges() const
void advanceCouplingScheme()
Advances the coupling schemes.
void computePartitions()
Communicate meshes and create partitions.
void setMeshEdge(std::string_view meshName, VertexID first, VertexID second)
Sets a mesh edge from vertex IDs.
std::vector< impl::PtrParticipant > _participants
Holds information about solvers participating in the coupled simulation.
bool requiresMeshConnectivityFor(std::string_view meshName) const
Checks if the given mesh requires connectivity.
void setMeshTriangles(std::string_view meshName, ::precice::span< const VertexID > vertices)
Sets multiple mesh triangles from vertex IDs.
int _executedReadMappings
Counts the amount of samples mapped in read mappings executed in the latest advance.
void clearStamplesOfChangedMeshes(MeshChanges totalMeshChanges)
Clears stample of changed meshes to make them consistent after the reinitialization.
void performDataActions(const std::set< action::Action::Timing > &timings)
Performs all data actions with given timing.
void setMeshTriangle(std::string_view meshName, VertexID first, VertexID second, VertexID third)
Sets mesh triangle from vertex IDs.
void configurePartitions(const m2n::M2NConfiguration::SharedPointer &m2nConfig)
Determines participants providing meshes to other participants.
double getMaxTimeStepSize() const
Get the maximum allowed time step size of the current window.
cplscheme::PtrCouplingScheme _couplingScheme
long int _numberAdvanceCalls
Counts calls to advance for plotting.
bool _allowsExperimental
Are experimental API calls allowed?
void handleDataBeforeAdvance(bool reachedTimeWindowEnd, double timeSteppedTo)
Completes everything data-related between adding time to and advancing the coupling scheme.
void setMeshTetrahedron(std::string_view meshName, VertexID first, VertexID second, VertexID third, VertexID fourth)
Set tetrahedron in 3D mesh from vertex ID.
void handleDataAfterAdvance(bool reachedTimeWindowEnd, bool isTimeWindowComplete, double timeSteppedTo, double timeAfterAdvance, const cplscheme::ImplicitData &receivedData)
Completes everything data-related after advancing the coupling scheme.
void samplizeWriteData(double time)
Creates a Stample at the given time for each write Data and zeros the buffers.
void setMeshQuads(std::string_view meshName, ::precice::span< const VertexID > vertices)
Sets multiple mesh quads from vertex IDs.
void getMeshVertexIDsAndCoordinates(std::string_view meshName, ::precice::span< VertexID > ids, ::precice::span< double > coordinates) const
getMeshVertexIDsAndCoordinates Iterates over the region of interest defined by bounding boxes and rea...
void closeCommunicationChannels(CloseChannels cc)
Syncs the primary ranks of all connected participants.
bool reinitHandshake(bool requestReinit) const
void trimSendDataAfter(double time)
Discards send (currently write) data of a participant after a given time when another iteration is re...
std::unique_ptr< profiling::Event > _solverInitEvent
ParticipantImpl(std::string_view participantName, std::string_view configurationFileName, int solverProcessIndex, int solverProcessSize, std::optional< void * > communicator)
Generic constructor for ParticipantImpl.
void advance(double computedTimeStepSize)
Advances preCICE after the solver has computed one time step.
void setMeshAccessRegion(std::string_view meshName, ::precice::span< const double > boundingBox) const
setMeshAccessRegion Define a region of interest on a received mesh (<receive-mesh ....
void writeData(std::string_view meshName, std::string_view dataName, ::precice::span< const VertexID > vertices, ::precice::span< const double > values)
Writes data to a mesh.
std::vector< profiling::Event > _userEvents
void handleExports(ExportTiming timing)
bool isCouplingOngoing() const
Checks if the coupled simulation is still ongoing.
State _state
The current State of the Participant.
void mapInitialWrittenData()
Computes, and performs write mappings of the initial data in initialize.
void setMeshVertices(std::string_view meshName, ::precice::span< const double > positions, ::precice::span< VertexID > ids)
Creates multiple mesh vertices.
void resetMesh(std::string_view meshName)
Removes all vertices and connectivity information from the mesh.
MappedSamples mappedSamples() const
Returns the amount of mapped read and write samples in the last call to advance.
void finalize()
Finalizes preCICE.
int _executedWriteMappings
Counts the amount of samples mapped in write mappings executed in the latest advance.
void configureM2Ns(const m2n::M2NConfiguration::SharedPointer &config)
void trimReadMappedData(double timeAfterAdvance, bool isTimeWindowComplete, const cplscheme::ImplicitData &fromData)
Removes samples in mapped to data connected to received data via a mapping.
std::map< std::string, m2n::BoundM2N > _m2ns
void reinitialize()
Reinitializes preCICE.
void setupWatcher()
Setup mesh watcher such as WatchPoints.
bool isTimeWindowComplete() const
Checks if the current coupling window is completed.
void setMeshEdges(std::string_view meshName, ::precice::span< const VertexID > vertices)
Sets multiple mesh edges from vertex IDs.
std::string _configHash
The hash of the configuration file used to configure this participant.
void setupCommunication()
Connect participants including repartitioning.
bool _allowsRemeshing
Are experimental remeshing API calls allowed?
VertexID setMeshVertex(std::string_view meshName, ::precice::span< const double > position)
Creates a mesh vertex.
void syncTimestep(double computedTimeStepSize)
Syncs the time step size between all ranks (all time steps sizes should be the same!...
void readData(std::string_view meshName, std::string_view dataName, ::precice::span< const VertexID > vertices, double relativeReadTime, ::precice::span< double > values) const
Reads data values from a mesh. Values correspond to a given point in time relative to the beginning o...
const mesh::Mesh & mesh(const std::string &meshName) const
Allows to access a registered mesh.
void compareBoundingBoxes()
Communicate bounding boxes and look for overlaps.
bool _waitInFinalize
Are participants waiting for each other in finalize?
void startProfilingSection(std::string_view eventName)
void initializeIntraCommunication()
Initializes intra-participant communication.
void resetWrittenData()
Resets written data.
void trimOldDataBefore(double time)
Discards data before the given time for all meshes and data known by this participant.
void mapAndReadData(std::string_view meshName, std::string_view dataName, ::precice::span< const double > coordinates, double relativeReadTime, ::precice::span< double > values) const
Reads data values from a mesh using a just-in-time data mapping. Values correspond to a given point i...
std::unique_ptr< profiling::Event > _solverAdvanceEvent
bool requiresUserDefinedAccessRegion(std::string_view meshName) const
void configure(std::string_view configurationFileName)
Configures the coupling interface from the given xml file.
void mapWrittenData(std::optional< double > after=std::nullopt)
Computes, and performs suitable write mappings either entirely or after given time.
void initialize()
Fully initializes preCICE and coupling data.
void writeAndMapData(std::string_view meshName, std::string_view dataName, ::precice::span< const double > coordinates, ::precice::span< const double > values)
Writes data values to a mesh using a just-in-time mapping (experimental).
void computeMappings(std::vector< MappingContext > &contexts, const std::string &mappingType)
Helper for mapWrittenData and mapReadData.
int getMeshDimensions(std::string_view meshName) const
Returns the spatial dimensionality of the given mesh.
Stores one Data object with related mesh. Context stores data to be read from and potentially provide...
bool hasSamples() const
Are there samples to read from?
void readValues(::precice::span< const VertexID > vertices, double time, ::precice::span< double > values) const
Samples data at a given point in time within the current time window for given indices.
void mapAndReadValues(::precice::span< const double > coordinates, double readTime, ::precice::span< double > values)
Forwards the just-in-time mapping API call for reading data to the data context.
Stores one Data object with related mesh. Context stores data to be written to and potentially provid...
void writeGradientsIntoDataBuffer(::precice::span< const VertexID > vertices, ::precice::span< const double > gradients)
Store gradients in _writeDataBuffer.
void writeAndMapValues(::precice::span< const double > coordinates, ::precice::span< const double > values)
Forwards the just-in-time mapping API call for writing data to the data context.
void writeValuesIntoDataBuffer(::precice::span< const VertexID > vertices, ::precice::span< const double > values)
Store values in _writeDataBuffer.
An M2N between participants with a configured direction.
Definition BoundM2N.hpp:12
std::string remoteName
Definition BoundM2N.hpp:31
std::string localName
Definition BoundM2N.hpp:30
Linear edge of a mesh, defined by two Vertex objects.
Definition Edge.hpp:15
Container and creator for meshes.
Definition Mesh.hpp:38
Triangle & createTriangle(Edge &edgeOne, Edge &edgeTwo, Edge &edgeThree)
Creates and initializes a Triangle object.
Definition Mesh.cpp:121
int getDimensions() const
Definition Mesh.cpp:100
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:220
std::size_t nVertices() const
Returns the number of vertices.
Definition Mesh.cpp:65
bool isValidVertexID(VertexID vertexID) const
Returns true if the given vertexID is valid.
Definition Mesh.cpp:230
Vertex & createVertex(const Eigen::Ref< const Eigen::VectorXd > &coords)
Creates and initializes a Vertex object.
Definition Mesh.cpp:105
void allocateDataValues()
Allocates memory for the vertex data values and corresponding gradient values.
Definition Mesh.cpp:235
Vertex of a mesh.
Definition Vertex.hpp:16
A partition that is provided by the participant.
A partition that is computed from a mesh received from another participant.
static EventRegistry & instance()
Returns the only instance (singleton) of the EventRegistry class.
void startBackend()
Create the file and starts the filestream if profiling is turned on.
void initialize(std::string_view applicationName, int rank=0, int size=1)
Sets the global start time.
void finalize()
Sets the global end time and flushes buffers.
void stop()
Stops a running event.
Definition Event.cpp:51
void addData(std::string_view key, int value)
Adds named integer data, associated to an event.
Definition Event.cpp:62
A C++ 11 implementation of the non-owning C++20 std::span type.
Definition span.hpp:284
constexpr pointer data() const noexcept
Definition span.hpp:500
PRECICE_SPAN_NODISCARD constexpr bool empty() const noexcept
Definition span.hpp:476
constexpr iterator begin() const noexcept
Definition span.hpp:503
constexpr iterator end() const noexcept
Definition span.hpp:505
constexpr size_type size() const noexcept
Definition span.hpp:469
static void barrier()
Synchronizes all ranks.
static void allreduceSum(precice::span< const double > sendData, precice::span< double > rcvData)
static bool isPrimary()
True if this process is running the primary rank.
Definition IntraComm.cpp:52
static void broadcast(bool &value)
static auto allSecondaryRanks()
Returns an iterable range over salve ranks [1, _size)
Definition IntraComm.hpp:37
static bool isParallel()
True if this process is running in parallel.
Definition IntraComm.cpp:62
static bool isSecondary()
True if this process is running a secondary rank.
Definition IntraComm.cpp:57
static com::PtrCommunication & getCommunication()
Intra-participant communication.
Definition IntraComm.hpp:31
static void configure(Rank rank, int size)
Configures the intra-participant communication.
Definition IntraComm.cpp:31
bool check(const K &name) const
Checks the status of a lock.
void add(Key name, bool state)
Adds a lock with a given state.
Definition MultiLock.hpp:41
void clear() noexcept
Removes all known locks.
Definition MultiLock.hpp:93
void lockAll() noexcept
Locks all known locks.
Definition MultiLock.hpp:62
bool checkAll() const noexcept
Checks whether all locks are locked.
void unlock(const K &name)
Unlocks a given lock.
Definition MultiLock.hpp:74
static void finalizeOrCleanupMPI()
Finalized a managed MPI environment or cleans up after an non-managed session.
Definition Parallel.cpp:230
static CommStatePtr current()
Returns an owning pointer to the current CommState.
Definition Parallel.cpp:147
static void initializeOrDetectMPI(std::optional< Communicator > userProvided=std::nullopt)
Definition Parallel.cpp:190
static void finalize()
Finalizes Petsc environment.
Definition Petsc.cpp:58
T current_path(T... args)
T distance(T... args)
T empty(T... args)
T end(T... args)
T find(T... args)
T get(T... args)
T make_unique(T... args)
std::string errorInvalidVertexID(int vid)
static constexpr auto errorInvalidVertexIDRange
void setMPIRank(int const rank)
void setParticipant(std::string const &participant)
ConvexityResult isConvexQuad(std::array< Eigen::VectorXd, 4 > coords)
Definition geometry.cpp:143
constexpr bool equals(const Eigen::MatrixBase< DerivedA > &A, const Eigen::MatrixBase< DerivedB > &B, double tolerance=NUMERICAL_ZERO_DIFFERENCE)
Compares two Eigen::MatrixBase for equality up to tolerance.
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type smallerEquals(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
constexpr double NUMERICAL_ZERO_DIFFERENCE
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greaterEquals(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
std::enable_if< std::is_arithmetic< Scalar >::value, bool >::type greater(Scalar A, Scalar B, Scalar tolerance=NUMERICAL_ZERO_DIFFERENCE)
std::array< Eigen::VectorXd, n > coordsFor(const Mesh &mesh, const std::array< int, n > &vertexIDs)
Given a mesh and an array of vertexIDS, this function returns an array of coordinates of the vertices...
Definition Utils.hpp:122
std::array< Vertex *, n > vertexPtrsFor(Mesh &mesh, const std::array< int, n > &vertexIDs)
Given a mesh and an array of vertexIDS, this function returns an array of pointers to vertices.
Definition Utils.hpp:111
std::size_t countVerticesInBoundingBox(mesh::PtrMesh mesh, const mesh::BoundingBox &bb)
Given a Mesh and a bounding box, counts all vertices within the bounding box.
Definition Utils.cpp:72
std::shared_ptr< Partition > PtrPartition
static constexpr FundamentalTag Fundamental
Convenience instance of the FundamentalTag.
Definition Event.hpp:19
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:21
auto reorder_array(const std::array< Index, n > &order, const std::array< T, n > &elements) -> std::array< T, n >
Reorders an array given an array of unique indices.
bool contained(const ELEMENT_T &element, const std::vector< ELEMENT_T > &vec)
Returns true, if given element is in vector, otherwise false.
Definition Helpers.hpp:38
std::pair< InputIt, InputIt > find_first_range(InputIt first, InputIt last, Predicate p)
Finds the first range in [first, last[ that fulfills a predicate.
bool unique_elements(const Container &c, BinaryPredicate p={})
Definition algorithm.hpp:63
auto make_array(Elements &&...elements) -> std::array< typename std::common_type< Elements... >::type, sizeof...(Elements)>
Function that generates an array from given elements.
Definition algorithm.hpp:50
std::string configure(XMLTag &tag, const precice::xml::ConfigurationContext &context, std::string_view configurationFilename)
Configures the given configuration from file configurationFilename.
Definition XMLTag.cpp:284
int VertexID
Definition Types.hpp:13
int Rank
Definition Types.hpp:37
T has_value(T... args)
T push_back(T... args)
T reset(T... args)
T size(T... args)
T sort(T... args)
T stable_partition(T... args)
Holds a data mapping and related information.
Stores a mesh and related objects and data.
mesh::PtrMesh mesh
Mesh holding the geometry data structure.
mapping::Mapping::MeshRequirement meshRequirement
Determines which mesh type has to be provided by the accessor.
void checkVerticesInsideAccessRegion(precice::span< const double > coordinates, const int meshDim, std::string_view functionName) const
std::shared_ptr< mesh::BoundingBox > userDefinedAccessRegion
std::vector< std::reference_wrapper< const mesh::Vertex > > filterVerticesToLocalAccessRegion(bool requiresBB) const
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21
T value(T... args)