preCICE
Loading...
Searching...
No Matches
ParticipantImpl.cpp
Go to the documentation of this file.
1#include <Eigen/Core>
2#include <algorithm>
3#include <array>
4#include <cmath>
5#include <deque>
6#include <filesystem>
7#include <functional>
8#include <iterator>
9#include <memory>
10#include <numeric>
11#include <optional>
12#include <ostream>
13#include <sstream>
14#include <string_view>
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"
66#include "precice/impl/versions.hpp"
67#include "profiling/Event.hpp"
71#include "utils/EigenIO.hpp"
72#include "utils/Helpers.hpp"
73#include "utils/IntraComm.hpp"
74#include "utils/Parallel.hpp"
75#include "utils/Petsc.hpp"
76#include "utils/algorithm.hpp"
77#include "utils/assertion.hpp"
78#include "xml/XMLTag.hpp"
79
81
82namespace precice::impl {
83
85 std::string_view participantName,
86 std::string_view configurationFileName,
87 int solverProcessIndex,
88 int solverProcessSize,
89 std::optional<void *> communicator)
90 : _accessorName(participantName),
91 _accessorProcessRank(solverProcessIndex),
92 _accessorCommunicatorSize(solverProcessSize)
93{
94
96 "This participant's name is an empty string. "
97 "When constructing a preCICE interface you need to pass the name of the "
98 "participant as first argument to the constructor.");
100
101 PRECICE_CHECK(!communicator || communicator.value() != nullptr,
102 "Passing \"nullptr\" as \"communicator\" to Participant constructor is not allowed. "
103 "Please use the Participant constructor without the \"communicator\" argument, if you don't want to pass an MPI communicator.");
105 "The solver process index needs to be a non-negative number, not: {}. "
106 "Please check the value given when constructing a preCICE interface.",
109 "The solver process size needs to be a positive number, not: {}. "
110 "Please check the value given when constructing a preCICE interface.",
113 "The solver process index, currently: {} needs to be smaller than the solver process size, currently: {}. "
114 "Please check the values given when constructing a preCICE interface.",
116
119 Event e("construction", profiling::Fundamental);
120
121 // Set the global communicator to the passed communicator.
122 // This is a noop if preCICE is not configured with MPI.
123#ifndef PRECICE_NO_MPI
124 Event e3("com.initializeMPI", profiling::Fundamental);
125 if (communicator.has_value()) {
126 auto commptr = static_cast<utils::Parallel::Communicator *>(communicator.value());
128 } else {
130 }
131
132 {
133 const auto currentRank = utils::Parallel::current()->rank();
134 PRECICE_CHECK(_accessorProcessRank == currentRank,
135 "The solver process index given in the preCICE interface constructor({}) does not match the rank of the passed MPI communicator ({}).",
136 _accessorProcessRank, currentRank);
137 const auto currentSize = utils::Parallel::current()->size();
139 "The solver process size given in the preCICE interface constructor({}) does not match the size of the passed MPI communicator ({}).",
140 _accessorCommunicatorSize, currentSize);
141 }
142 e3.stop();
143#else
144 PRECICE_WARN_IF(communicator.has_value(), "preCICE was configured without MPI but you passed an MPI communicator. preCICE ignores the communicator and continues.");
145#endif
146
147 Event e1("configure", profiling::Fundamental);
148 configure(configurationFileName);
149 e1.stop();
150
151 // Backend settings have been configured
152 Event e2("startProfilingBackend");
154 e2.stop();
155
156 PRECICE_DEBUG("Initialize intra-participant communication");
159 }
160
161 e.stop();
162 _solverInitEvent = std::make_unique<profiling::Event>("solver.initialize", profiling::Fundamental, profiling::Synchronize);
163}
164
166{
167 if (_state != State::Finalized) {
168 PRECICE_INFO("Implicitly finalizing in destructor");
169 finalize();
170 }
171}
172
174 std::string_view configurationFileName)
175{
176
183 _configHash = xml::configure(config.getXMLTag(), context, configurationFileName);
184 if (_accessorProcessRank == 0) {
185 PRECICE_INFO("This is preCICE version {}", PRECICE_VERSION);
186 PRECICE_INFO("Revision info: {}", precice::preciceRevision);
187 constexpr std::string_view buildTypeStr = "Build type: "
188#ifndef NDEBUG
189 "Debug"
190#else // NDEBUG
191 "Release"
192#ifndef PRECICE_NO_DEBUG_LOG
193 " + debug log"
194#else
195 " (without debug log)"
196#endif
197#ifndef PRECICE_NO_TRACE_LOG
198 " + trace log"
199#endif
200#ifndef PRECICE_NO_ASSERTIONS
201 " + assertions"
202#endif
203#endif // NDEBUG
204 ;
205 PRECICE_INFO(buildTypeStr);
206 try {
207 PRECICE_INFO("Working directory \"{}\"", std::filesystem::current_path().string());
208 } catch (std::filesystem::filesystem_error &fse) {
209 PRECICE_INFO("Working directory unknown due to error \"{}\"", fse.what());
210 }
211 PRECICE_INFO("Configuring preCICE with configuration \"{}\"", configurationFileName);
212 PRECICE_INFO("I am participant \"{}\"", _accessorName);
213 }
214
216
217 PRECICE_CHECK(config.getParticipantConfiguration()->nParticipants() > 1,
218 "In the preCICE configuration, only one participant is defined. "
219 "One participant makes no coupled simulation. "
220 "Please add at least another one.");
221
222 _allowsExperimental = config.allowsExperimental();
223 _allowsRemeshing = config.allowsRemeshing();
224 _waitInFinalize = config.waitInFinalize();
226 _participants = config.getParticipantConfiguration()->getParticipants();
227 _m2ns = config.getBoundM2NsFor(_accessorName);
228 config.configurePartitionsFor(_accessorName);
229 _couplingScheme = config.getCouplingSchemeConfiguration()->getCouplingScheme(_accessorName);
230
231 PRECICE_ASSERT(_accessorCommunicatorSize == 1 || _accessor->useIntraComm(),
232 "A parallel participant needs an intra-participant communication");
233 PRECICE_CHECK(not(_accessorCommunicatorSize == 1 && _accessor->useIntraComm()),
234 "You cannot use an intra-participant communication with a serial participant. "
235 "If you do not know exactly what an intra-participant communication is and why you want to use it "
236 "you probably just want to remove the intraComm tag from the preCICE configuration.");
237
239
240 // Register all MeshIds to the lock, but unlock them straight away as
241 // writing is allowed after configuration.
242 _meshLock.clear();
243 for (const auto &variant : _accessor->usedMeshContexts()) {
244 _meshLock.add(getMesh(variant).getName(), false);
245 }
246}
247
249{
251 PRECICE_CHECK(_state != State::Finalized, "initialize() cannot be called after finalize().");
252 PRECICE_CHECK(_state != State::Initialized, "initialize() may only be called once.");
253 PRECICE_ASSERT(not _couplingScheme->isInitialized());
254
256 PRECICE_CHECK(not failedToInitialize,
257 "Initial data has to be written to preCICE before calling initialize(). "
258 "After defining your mesh, call requiresInitialData() to check if the participant is required to write initial data using the writeData() function.");
259
260 // Enforce that all user-created events are stopped to prevent incorrect nesting.
261 PRECICE_CHECK(_userEvents.empty(), "There are unstopped user defined events. Please stop them using stopLastProfilingSection() before calling initialize().");
262
263 _solverInitEvent.reset();
265
266 for (const auto &context : _accessor->providedMeshContexts()) {
267 e.addData("meshSize" + context.mesh->getName(), context.mesh->nVertices());
268 }
269
271 setupWatcher();
272
273 _meshLock.lockAll();
274
275 for (auto &context : _accessor->writeDataContexts()) {
276 const double startTime = 0.0;
277 context.storeBufferedData(startTime);
278 }
279
282
283 PRECICE_DEBUG("Initialize coupling schemes");
284 Event e1("initalizeCouplingScheme", profiling::Fundamental);
285 _couplingScheme->initialize();
286 e1.stop();
287
290
292
294
295 e.stop();
296
298 PRECICE_INFO(_couplingScheme->printCouplingState());
299 _solverAdvanceEvent = std::make_unique<profiling::Event>("solver.advance", profiling::Fundamental, profiling::Synchronize);
300}
301
303{
306 PRECICE_INFO("Reinitializing Participant");
307 Event e("reinitialize", profiling::Fundamental);
309
310 for (const auto &context : _accessor->providedMeshContexts()) {
311 e.addData("meshSize" + context.mesh->getName(), context.mesh->nVertices());
312 }
313
315 setupWatcher();
316
317 PRECICE_DEBUG("Reinitialize coupling schemes");
318 _couplingScheme->reinitialize();
319}
320
322{
324
325 // TODO only preprocess changed meshes
326 PRECICE_DEBUG("Preprocessing provided meshes");
327 for (auto &context : _accessor->providedMeshContexts()) {
328 auto &mesh = *context.mesh;
329 Event e("preprocess." + mesh.getName());
330 mesh.preprocess();
331 }
332
333 // Setup communication
334
335 PRECICE_INFO("Setting up primary communication to coupling partner/s");
336 Event e2("connectPrimaries");
337 for (auto &m2nPair : _m2ns) {
338 auto &bm2n = m2nPair.second;
339 bool requesting = bm2n.isRequesting;
340 if (bm2n.m2n->isConnected()) {
341 PRECICE_DEBUG("Primary connection {} {} already connected.", (requesting ? "from" : "to"), bm2n.remoteName);
342 } else {
343 PRECICE_DEBUG((requesting ? "Awaiting primary connection from {}" : "Establishing primary connection to {}"), bm2n.remoteName);
344 bm2n.prepareEstablishment();
345 bm2n.connectPrimaryRanks(_configHash);
346 PRECICE_DEBUG("Established primary connection {} {}", (requesting ? "from " : "to "), bm2n.remoteName);
347 }
348 }
349 e2.stop();
350
351 PRECICE_INFO("Primary ranks are connected");
352
353 Event e3("repartitioning");
354 // clears the mappings as well (see clearMappings)
356
357 PRECICE_INFO("Setting up preliminary secondary communication to coupling partner/s");
358 for (auto &m2nPair : _m2ns) {
359 auto &bm2n = m2nPair.second;
360 bm2n.preConnectSecondaryRanks();
361 }
362
364 e3.stop();
365
366 PRECICE_INFO("Setting up secondary communication to coupling partner/s");
367 Event e4("connectSecondaries");
368 for (auto &m2nPair : _m2ns) {
369 auto &bm2n = m2nPair.second;
370 bm2n.connectSecondaryRanks();
371 PRECICE_DEBUG("Established secondary connection {} {}", (bm2n.isRequesting ? "from " : "to "), bm2n.remoteName);
372 }
373 PRECICE_INFO("Secondary ranks are connected");
374
375 for (auto &m2nPair : _m2ns) {
376 m2nPair.second.cleanupEstablishment();
377 }
378}
379
381{
383 PRECICE_DEBUG("Initialize watchpoints");
384 for (PtrWatchPoint &watchPoint : _accessor->watchPoints()) {
385 watchPoint->initialize();
386 }
387 for (PtrWatchIntegral &watchIntegral : _accessor->watchIntegrals()) {
388 watchIntegral->initialize();
389 }
390}
391
393 double computedTimeStepSize)
394{
395
396 PRECICE_TRACE(computedTimeStepSize);
397
398 // Enforce that all user-created events are stopped to prevent incorrect nesting.
399 PRECICE_CHECK(_userEvents.empty(), "There are unstopped user defined events. Please stop them using stopLastProfilingSection() before calling advance().");
400
401 // Events for the solver time, stopped when we enter, restarted when we leave advance
402 PRECICE_ASSERT(_solverAdvanceEvent, "The advance event is created in initialize");
403 _solverAdvanceEvent->stop();
404
406
407 PRECICE_CHECK(_state != State::Constructed, "initialize() has to be called before advance().");
408 PRECICE_CHECK(_state != State::Finalized, "advance() cannot be called after finalize().");
409 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before advance().");
410 PRECICE_ASSERT(_couplingScheme->isInitialized());
411 PRECICE_CHECK(isCouplingOngoing(), "advance() cannot be called when isCouplingOngoing() returns false.");
412
413 // validating computed time step
414 PRECICE_CHECK(std::isfinite(computedTimeStepSize), "advance() cannot be called with an infinite time step size.");
415 PRECICE_CHECK(!math::equals(computedTimeStepSize, 0.0), "advance() cannot be called with a time step size of 0.");
416 PRECICE_CHECK(computedTimeStepSize > 0.0, "advance() cannot be called with a negative time step size {}.", computedTimeStepSize);
417
419
420#ifndef NDEBUG
421 PRECICE_DEBUG("Synchronize time step size");
423 syncTimestep(computedTimeStepSize);
424 }
425#endif
426
427 // Update the coupling scheme time state. Necessary to get correct remainder.
428 const bool isAtWindowEnd = _couplingScheme->addComputedTime(computedTimeStepSize);
429
430 if (_allowsRemeshing) {
431 if (isAtWindowEnd) {
432 auto totalMeshChanges = getTotalMeshChanges();
433 clearStamplesOfChangedMeshes(totalMeshChanges);
434
435 int sumOfChanges = std::accumulate(totalMeshChanges.begin(), totalMeshChanges.end(), 0);
436 if (reinitHandshake(sumOfChanges)) {
437 reinitialize();
438 }
439 } else {
440 PRECICE_CHECK(_meshLock.checkAll(), "The time window needs to end after remeshing.");
441 }
442 }
443
444 const double timeSteppedTo = _couplingScheme->getTime();
445 const auto dataToReceive = _couplingScheme->implicitDataToReceive();
446
447 handleDataBeforeAdvance(isAtWindowEnd, timeSteppedTo);
448
450
451 // In clase if an implicit scheme, this may be before timeSteppedTo
452 const double timeAfterAdvance = _couplingScheme->getTime();
453 const bool timeWindowComplete = _couplingScheme->isTimeWindowComplete();
454
455 handleDataAfterAdvance(isAtWindowEnd, timeWindowComplete, timeSteppedTo, timeAfterAdvance, dataToReceive);
456
457 PRECICE_INFO(_couplingScheme->printCouplingState());
458
459 PRECICE_DEBUG("Mapped {} samples in write mappings and {} samples in read mappings",
461
462 _meshLock.lockAll();
463
464 e.stop();
465 _solverAdvanceEvent->start();
466}
467
468void ParticipantImpl::handleDataBeforeAdvance(bool reachedTimeWindowEnd, double timeSteppedTo)
469{
470 // We only have to care about write data, in case substeps are enabled
471 // OR we are at the end of a timewindow, otherwise, we simply erase
472 // them as they have no relevance for the coupling (without time
473 // interpolation, only the time window end is relevant), the resetting
474 // happens regardless of the if-condition.
475 if (reachedTimeWindowEnd || _couplingScheme->requiresSubsteps()) {
476
477 // Here, we add the written data to the waveform storage. In the
478 // mapWrittenData, we then take samples from the storage and execute
479 // the mapping using waveform samples on the (for write mappings) "to"
480 // side.
481 samplizeWriteData(timeSteppedTo);
482 }
483
485
486 // Reset mapping counters here to cover subcycling
489
490 if (reachedTimeWindowEnd) {
491 mapWrittenData(_couplingScheme->getTimeWindowStart());
493 }
494}
495
496void ParticipantImpl::handleDataAfterAdvance(bool reachedTimeWindowEnd, bool isTimeWindowComplete, double timeSteppedTo, double timeAfterAdvance, const cplscheme::ImplicitData &receivedData)
497{
498 if (!reachedTimeWindowEnd) {
499 // We are subcycling
500 return;
501 }
502
504 // Move to next time window
505 PRECICE_ASSERT(math::greaterEquals(timeAfterAdvance, timeSteppedTo), "We must have stayed or moved forwards in time (min-time-step-size).", timeAfterAdvance, timeSteppedTo);
506
507 // As we move forward, there may now be old samples lying around
508 trimOldDataBefore(_couplingScheme->getTimeWindowStart());
509 } else {
510 // We are iterating
511 PRECICE_ASSERT(math::greater(timeSteppedTo, timeAfterAdvance), "We must have moved back in time!");
512
513 trimSendDataAfter(timeAfterAdvance);
514 }
515
516 if (reachedTimeWindowEnd) {
517 trimReadMappedData(timeAfterAdvance, isTimeWindowComplete, receivedData);
518 mapReadData();
520 }
521
522 // Required for implicit coupling
523 for (auto &context : _accessor->readDataContexts()) {
524 context.invalidateMappingCache();
525 }
526
527 // Strictly speaking, the write direction is not relevant here, but we will add it for the sake of completenss
528 for (auto &context : _accessor->writeDataContexts()) {
529 context.invalidateMappingCache();
530 }
531
533 // Reset initial guesses for iterative mappings
534 for (auto &context : _accessor->readDataContexts()) {
535 context.resetInitialGuesses();
536 }
537 for (auto &context : _accessor->writeDataContexts()) {
538 context.resetInitialGuesses();
539 }
540 }
541
543}
544
546{
547 // store buffered write data in sample storage and reset the buffer
548 for (auto &context : _accessor->writeDataContexts()) {
549
550 // Finalize conservative write mapping, later we reset
551 // the buffer in resetWrittenData
552
553 // Note that "samplizeWriteData" operates on _providedData of the
554 // DataContext, which is for just-in-time mappings the data we write
555 // on the received mesh.
556 // For just-in-time mappings, the _providedData should contain by now
557 // the "just-in-time" mapped data. However, it would be wasteful to
558 // execute expensive parts (in particular solving the RBF systems)
559 // for each writeAndMapData call. Thus, we create a DataCache during
560 // the writeAndMapData API calls, which contains pre-processed data
561 // values. Here, we now need to finalize the just-in-time mappings,
562 // before we can add it to the waveform buffer.
563 // For now, this only applies to just-in-time write mappings
564
565 context.completeJustInTimeMapping();
566 context.storeBufferedData(time);
567 }
568}
569
571{
572 for (auto &variant : _accessor->usedMeshContexts()) {
573 auto &mesh = getMesh(variant);
574 for (const auto &name : mesh.availableData()) {
575 mesh.data(name)->waveform().trimBefore(time);
576 }
577 }
578}
579
581{
582 for (auto &context : _accessor->writeDataContexts()) {
583 context.trimAfter(time);
584 }
585}
586
588{
590 PRECICE_CHECK(_state != State::Finalized, "finalize() may only be called once.");
591
592 // First we gracefully stop all existing user events and finally the last solver.advance event
593 while (!_userEvents.empty()) {
594 // Ensure reverse destruction order for correct nesting
595 _userEvents.pop_back();
596 }
597 _solverAdvanceEvent.reset();
598
599 Event e("finalize", profiling::Fundamental);
600
601 if (_state == State::Initialized) {
602
603 PRECICE_ASSERT(_couplingScheme->isInitialized());
604 PRECICE_DEBUG("Finalize coupling scheme");
605 _couplingScheme->finalize();
606
608 }
609
610 // Release ownership
611 _couplingScheme.reset();
612 _participants.clear();
613 _accessor.reset();
614
615 // Close Connections
616 PRECICE_DEBUG("Close intra-participant communication");
618 utils::IntraComm::getCommunication()->closeConnection();
620 }
621 _m2ns.clear();
622
623 // Stop and print Event logging
624 e.stop();
625
626 // Finalize PETSc and Events first
628// This will lead to issues if we call finalize afterwards again
629#if !defined(PRECICE_NO_GINKGO) || !defined(PRECICE_NO_KOKKOS_KERNELS)
631#endif
633
634 // Finally clear events and finalize MPI
637}
638
639int ParticipantImpl::getMeshDimensions(std::string_view meshName) const
640{
641 PRECICE_TRACE(meshName);
643 return _accessor->meshContext(meshName).mesh->getDimensions();
644}
645
646int ParticipantImpl::getDataDimensions(std::string_view meshName, std::string_view dataName) const
647{
648 PRECICE_TRACE(meshName, dataName);
650 PRECICE_VALIDATE_DATA_NAME(meshName, dataName);
651 return _accessor->meshContext(meshName).mesh->data(dataName)->getDimensions();
652}
653
655{
657 PRECICE_CHECK(_state != State::Finalized, "isCouplingOngoing() cannot be called after finalize().");
658 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before isCouplingOngoing() can be evaluated.");
659 return _couplingScheme->isCouplingOngoing();
660}
661
663{
665 PRECICE_CHECK(_state != State::Constructed, "initialize() has to be called before isTimeWindowComplete().");
666 PRECICE_CHECK(_state != State::Finalized, "isTimeWindowComplete() cannot be called after finalize().");
667 return _couplingScheme->isTimeWindowComplete();
668}
669
671{
672 PRECICE_CHECK(_state != State::Finalized, "getMaxTimeStepSize() cannot be called after finalize().");
673 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before getMaxTimeStepSize() can be evaluated.");
674 const double nextTimeStepSize = _couplingScheme->getNextTimeStepMaxSize();
675 // PRECICE_ASSERT(!math::equals(nextTimeStepSize, 0.0), nextTimeStepSize); // @todo requires https://github.com/precice/precice/issues/1904
676 // PRECICE_ASSERT(math::greater(nextTimeStepSize, 0.0), nextTimeStepSize); // @todo requires https://github.com/precice/precice/issues/1904
677
678 // safeguard needed because _couplingScheme->getNextTimeStepMaxSize() returns 0, if not isCouplingOngoing()
679 // actual case where we want to warn the user
681 isCouplingOngoing() && not math::greater(nextTimeStepSize, 0.0, 100 * math::NUMERICAL_ZERO_DIFFERENCE),
682 "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.",
683 nextTimeStepSize);
684 return nextTimeStepSize;
685}
686
688{
690 PRECICE_CHECK(_state == State::Constructed, "requiresInitialData() has to be called before initialize().");
692 if (required) {
694 }
695 return required;
696}
697
699{
701 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before requiresWritingCheckpoint().");
703 if (required) {
705 }
706 return required;
707}
708
710{
712 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before requiresReadingCheckpoint().");
714 if (required) {
716 }
717 return required;
718}
719
720bool ParticipantImpl::requiresMeshConnectivityFor(std::string_view meshName) const
721{
723 MeshContext &context = _accessor->meshContext(meshName);
725}
726
727bool ParticipantImpl::requiresGradientDataFor(std::string_view meshName,
728 std::string_view dataName) const
729{
730 PRECICE_VALIDATE_DATA_NAME(meshName, dataName);
731 // Read data never requires gradients
732 if (!_accessor->isDataWrite(meshName, dataName))
733 return false;
734
735 WriteDataContext &context = _accessor->writeDataContext(meshName, dataName);
736 return context.hasGradient();
737}
738
740 std::string_view meshName) const
741{
742 PRECICE_TRACE(meshName);
743 PRECICE_REQUIRE_MESH_USE(meshName);
744 // In case we access received mesh data: check, if the requested mesh data has already been received.
745 // Otherwise, the function call doesn't make any sense
746 PRECICE_CHECK((_state == State::Initialized) || _accessor->isMeshProvided(meshName),
747 "initialize() has to be called before accessing data of the received mesh \"{}\" on participant \"{}\".",
748 meshName, _accessor->getName());
749
750 // Returns true if we have api access configured and we run in parallel and have a received mesh
751 if (_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName)) {
752 auto &receivedContext = _accessor->receivedMeshContext(meshName);
753 if (receivedContext.userDefinedAccessRegion || requiresUserDefinedAccessRegion(meshName)) {
754 // filter nVertices to the actual number of vertices queried by the user
755 PRECICE_CHECK(receivedContext.userDefinedAccessRegion, "The function getMeshVertexSize was called on the received mesh \"{0}\", "
756 "but no access region was defined although this is necessary for parallel runs. "
757 "Please define an access region using \"setMeshAccessRegion()\" before calling \"getMeshVertexSize()\".",
758 meshName);
759
760 auto result = mesh::countVerticesInBoundingBox(receivedContext.mesh, *receivedContext.userDefinedAccessRegion);
761
762 PRECICE_DEBUG("Filtered {} of {} vertices out on mesh {} due to the local access region. Mesh size in the access region: {}", receivedContext.mesh->nVertices() - result, receivedContext.mesh->nVertices(), meshName, result);
763 return result;
764 }
765 }
766 // For provided meshes and in case the api-access was not configured, we return here all vertices
767 PRECICE_WARN_IF(_accessor->isMeshReceived(meshName) && !_accessor->isDirectAccessAllowed(meshName),
768 "You are calling \"getMeshVertexSize()\" on a received mesh without api-access enabled (<receive-mesh name=\"{0}\" ... api-access=\"false\"/>). "
769 "Note that enabling api-access is required for this function to work properly with direct mesh access and just-in-time mappings.",
770 meshName);
771 return _accessor->meshContext(meshName).mesh->nVertices();
772}
773
776 std::string_view meshName)
777{
779 PRECICE_CHECK(_allowsRemeshing, "Cannot reset meshes. This feature needs to be enabled using <precice-configuration experimental=\"1\" allow-remeshing=\"1\">.");
780 PRECICE_CHECK(_state == State::Initialized, "initialize() has to be called before resetMesh().");
781 PRECICE_TRACE(meshName);
783 PRECICE_CHECK(_couplingScheme->isCouplingOngoing(), "Cannot remesh after the last time window has been completed.");
784 PRECICE_CHECK(_couplingScheme->isTimeWindowComplete(), "Cannot remesh while subcycling or iterating. Remeshing is only allowed when the time window is completed.");
785 impl::MeshContext &context = _accessor->meshContext(meshName);
786
787 PRECICE_DEBUG("Clear mesh positions for mesh \"{}\"", context.mesh->getName());
788 _meshLock.unlock(meshName);
789 context.mesh->clear();
790}
791
793 std::string_view meshName,
795{
796 PRECICE_TRACE(meshName);
798 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
799 auto &mesh = *context.mesh;
800 PRECICE_CHECK(position.size() == static_cast<unsigned long>(mesh.getDimensions()),
801 "Cannot set vertex for mesh \"{}\". Expected {} position components but found {}.", meshName, mesh.getDimensions(), position.size());
802 Event e{fmt::format("setMeshVertex.{}", meshName), profiling::API};
803 auto index = mesh.createVertex(Eigen::Map<const Eigen::VectorXd>{position.data(), mesh.getDimensions()}).getID();
804 mesh.allocateDataValues();
805
806 const auto newSize = mesh.nVertices();
807 for (auto &context : _accessor->writeDataContexts()) {
808 if (context.getMeshName() == mesh.getName()) {
809 context.resizeBufferTo(newSize);
810 }
811 }
812
813 return index;
814}
815
817 std::string_view meshName,
820{
821 PRECICE_TRACE(meshName, positions.size(), ids.size());
823 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
824 auto &mesh = *context.mesh;
825
826 const auto meshDims = mesh.getDimensions();
827 const auto expectedPositionSize = ids.size() * meshDims;
828 PRECICE_CHECK(positions.size() == expectedPositionSize,
829 "Input sizes are inconsistent attempting to set vertices on {}D mesh \"{}\". "
830 "You passed {} vertex indices and {} position components, but we expected {} position components ({} x {}).",
831 meshDims, meshName, ids.size(), positions.size(), expectedPositionSize, ids.size(), meshDims);
832
833 Event e{fmt::format("setMeshVertices.{}", meshName), profiling::API};
834 const Eigen::Map<const Eigen::MatrixXd> posMatrix{
835 positions.data(), mesh.getDimensions(), static_cast<EIGEN_DEFAULT_DENSE_INDEX_TYPE>(ids.size())};
836 for (unsigned long i = 0; i < ids.size(); ++i) {
837 ids[i] = mesh.createVertex(posMatrix.col(i)).getID();
838 }
839 mesh.allocateDataValues();
840
841 const auto newSize = mesh.nVertices();
842 for (auto &context : _accessor->writeDataContexts()) {
843 if (context.getMeshName() == mesh.getName()) {
844 context.resizeBufferTo(newSize);
845 }
846 }
847}
848
850 std::string_view meshName,
851 VertexID first,
852 VertexID second)
853{
854 PRECICE_TRACE(meshName, first, second);
856 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
858 return;
859 }
860
861 mesh::Mesh &mesh = *context.mesh;
863 PRECICE_CHECK(mesh.isValidVertexID(first), errorInvalidVertexID(first));
864 PRECICE_CHECK(mesh.isValidVertexID(second), errorInvalidVertexID(second));
865 Event e{fmt::format("setMeshEdge.{}", meshName), profiling::API};
866 mesh::Vertex &v0 = mesh.vertex(first);
867 mesh::Vertex &v1 = mesh.vertex(second);
868 mesh.createEdge(v0, v1);
869}
870
872 std::string_view meshName,
874{
875 PRECICE_TRACE(meshName, vertices.size());
877 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
879 return;
880 }
881
882 mesh::Mesh &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::API};
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, second, third);
915 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
917 return;
918 }
919
920 mesh::Mesh &mesh = *context.mesh;
922 PRECICE_CHECK(mesh.isValidVertexID(first), errorInvalidVertexID(first));
923 PRECICE_CHECK(mesh.isValidVertexID(second), errorInvalidVertexID(second));
924 PRECICE_CHECK(mesh.isValidVertexID(third), errorInvalidVertexID(third));
926 "setMeshTriangle() was called with repeated Vertex IDs ({}, {}, {}).",
927 first, second, third);
928
929 mesh::Vertex &A = mesh.vertex(first);
930 mesh::Vertex &B = mesh.vertex(second);
931 mesh::Vertex &C = mesh.vertex(third);
932
933 mesh.createTriangle(A, B, C);
934}
935
937 std::string_view meshName,
939{
940 PRECICE_TRACE(meshName, vertices.size());
942 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
944 return;
945 }
946
947 mesh::Mesh &mesh = *context.mesh;
948 PRECICE_CHECK(vertices.size() % 3 == 0,
949 "Cannot interpret passed vertex IDs attempting to set triangles of mesh \"{}\" . "
950 "You passed {} vertex indices, which isn't dividable by 3.",
951 meshName, vertices.size());
952 {
953 auto end = vertices.end();
954 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
955 return !mesh.isValidVertexID(vid);
956 });
957 PRECICE_CHECK(first == end,
959 std::distance(vertices.begin(), first),
960 std::distance(vertices.begin(), last));
961 }
962
963 Event e{fmt::format("setMeshTriangles.{}", meshName), profiling::API};
964
965 for (unsigned long i = 0; i < vertices.size() / 3; ++i) {
966 auto aid = vertices[3 * i];
967 auto bid = vertices[3 * i + 1];
968 auto cid = vertices[3 * i + 2];
969 mesh.createTriangle(mesh.vertex(aid),
970 mesh.vertex(bid),
971 mesh.vertex(cid));
972 }
973}
974
976 std::string_view meshName,
977 VertexID first,
978 VertexID second,
979 VertexID third,
980 VertexID fourth)
981{
982 PRECICE_TRACE(meshName, first,
983 second, third, fourth);
985 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
987 return;
988 }
989
990 PRECICE_ASSERT(context.mesh);
991 mesh::Mesh &mesh = *context.mesh;
993 PRECICE_CHECK(mesh.isValidVertexID(first), errorInvalidVertexID(first));
994 PRECICE_CHECK(mesh.isValidVertexID(second), errorInvalidVertexID(second));
995 PRECICE_CHECK(mesh.isValidVertexID(third), errorInvalidVertexID(third));
996 PRECICE_CHECK(mesh.isValidVertexID(fourth), errorInvalidVertexID(fourth));
997
998 auto vertexIDs = utils::make_array(first, second, third, fourth);
999 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.");
1000
1001 auto coords = mesh::coordsFor(mesh, vertexIDs);
1003 "The four vertices that form the quad are not unique. The resulting shape may be a point, line or triangle. "
1004 "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.");
1005
1006 auto convexity = math::geometry::isConvexQuad(coords);
1007 PRECICE_CHECK(convexity.convex, "The given quad is not convex. "
1008 "Please check that the adapter send the four correct vertices or that the interface is composed of quads.");
1009 auto reordered = utils::reorder_array(convexity.vertexOrder, mesh::vertexPtrsFor(mesh, vertexIDs));
1010
1011 Event e{fmt::format("setMeshQuad.{}", meshName), profiling::API};
1012
1013 // Vertices are now in the order: V0-V1-V2-V3-V0.
1014 // Use the shortest diagonal to split the quad into 2 triangles.
1015 // Vertices are now in V0-V1-V2-V3-V0 order. The new edge, e[4] is either 0-2 or 1-3
1016 double distance02 = (reordered[0]->getCoords() - reordered[2]->getCoords()).norm();
1017 double distance13 = (reordered[1]->getCoords() - reordered[3]->getCoords()).norm();
1018
1019 // The new edge, e[4], is the shortest diagonal of the quad
1020 if (distance02 <= distance13) {
1021 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[1]);
1022 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[3]);
1023 } else {
1024 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[0]);
1025 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[2]);
1026 }
1027}
1028
1030 std::string_view meshName,
1032{
1033 PRECICE_TRACE(meshName, vertices.size());
1035 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
1037 return;
1038 }
1039
1040 mesh::Mesh &mesh = *context.mesh;
1041 PRECICE_CHECK(vertices.size() % 4 == 0,
1042 "Cannot interpret passed vertex IDs attempting to set quads of mesh \"{}\" . "
1043 "You passed {} vertex indices, which isn't dividable by 4.",
1044 meshName, vertices.size());
1045 {
1046 auto end = vertices.end();
1047 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
1048 return !mesh.isValidVertexID(vid);
1049 });
1050 PRECICE_CHECK(first == end,
1052 std::distance(vertices.begin(), first),
1053 std::distance(vertices.begin(), last));
1054 }
1055
1056 for (unsigned long i = 0; i < vertices.size() / 4; ++i) {
1057 auto aid = vertices[4 * i];
1058 auto bid = vertices[4 * i + 1];
1059 auto cid = vertices[4 * i + 2];
1060 auto did = vertices[4 * i + 3];
1061
1062 auto vertexIDs = utils::make_array(aid, bid, cid, did);
1063 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);
1064
1065 auto coords = mesh::coordsFor(mesh, vertexIDs);
1067 "The four vertices that form the quad nr {} are not unique. The resulting shape may be a point, line or triangle. "
1068 "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.",
1069 i);
1070
1071 auto convexity = math::geometry::isConvexQuad(coords);
1072 PRECICE_CHECK(convexity.convex, "The given quad nr {} is not convex. "
1073 "Please check that the adapter send the four correct vertices or that the interface is composed of quads.",
1074 i);
1075 auto reordered = utils::reorder_array(convexity.vertexOrder, mesh::vertexPtrsFor(mesh, vertexIDs));
1076
1077 Event e{fmt::format("setMeshQuads.{}", meshName), profiling::API};
1078
1079 // Use the shortest diagonal to split the quad into 2 triangles.
1080 // Vertices are now in V0-V1-V2-V3-V0 order. The new edge, e[4] is either 0-2 or 1-3
1081 double distance02 = (reordered[0]->getCoords() - reordered[2]->getCoords()).norm();
1082 double distance13 = (reordered[1]->getCoords() - reordered[3]->getCoords()).norm();
1083
1084 if (distance02 <= distance13) {
1085 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[1]);
1086 mesh.createTriangle(*reordered[0], *reordered[2], *reordered[3]);
1087 } else {
1088 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[0]);
1089 mesh.createTriangle(*reordered[1], *reordered[3], *reordered[2]);
1090 }
1091 }
1092}
1093
1095 std::string_view meshName,
1096 VertexID first,
1097 VertexID second,
1098 VertexID third,
1099 VertexID fourth)
1100{
1101 PRECICE_TRACE(meshName, first, second, third, fourth);
1103 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
1104 PRECICE_CHECK(context.mesh->getDimensions() == 3, "setMeshTetrahedron is only possible for 3D meshes. "
1105 "Please set the mesh dimension to 3 in the preCICE configuration file.");
1107 return;
1108 }
1109
1110 Event e{fmt::format("setMeshTetrahedron.{}", meshName), profiling::API};
1111
1112 mesh::Mesh &mesh = *context.mesh;
1114 PRECICE_CHECK(mesh.isValidVertexID(first), errorInvalidVertexID(first));
1115 PRECICE_CHECK(mesh.isValidVertexID(second), errorInvalidVertexID(second));
1116 PRECICE_CHECK(mesh.isValidVertexID(third), errorInvalidVertexID(third));
1117 PRECICE_CHECK(mesh.isValidVertexID(fourth), errorInvalidVertexID(fourth));
1118 mesh::Vertex &A = mesh.vertex(first);
1119 mesh::Vertex &B = mesh.vertex(second);
1120 mesh::Vertex &C = mesh.vertex(third);
1121 mesh::Vertex &D = mesh.vertex(fourth);
1122
1123 mesh.createTetrahedron(A, B, C, D);
1124}
1125
1127 std::string_view meshName,
1129{
1130 PRECICE_TRACE(meshName, vertices.size());
1132 ProvidedMeshContext &context = _accessor->providedMeshContext(meshName);
1133 PRECICE_CHECK(context.mesh->getDimensions() == 3, "setMeshTetrahedron is only possible for 3D meshes. "
1134 "Please set the mesh dimension to 3 in the preCICE configuration file.");
1136 return;
1137 }
1138
1139 mesh::Mesh &mesh = *context.mesh;
1140 PRECICE_CHECK(vertices.size() % 4 == 0,
1141 "Cannot interpret passed vertex IDs attempting to set quads of mesh \"{}\" . "
1142 "You passed {} vertex indices, which isn't dividable by 4.",
1143 meshName, vertices.size());
1144 {
1145 auto end = vertices.end();
1146 auto [first, last] = utils::find_first_range(vertices.begin(), end, [&mesh](VertexID vid) {
1147 return !mesh.isValidVertexID(vid);
1148 });
1149 PRECICE_CHECK(first == end,
1151 std::distance(vertices.begin(), first),
1152 std::distance(vertices.begin(), last));
1153 }
1154
1155 Event e{fmt::format("setMeshTetrahedra.{}", meshName), profiling::API};
1156
1157 for (unsigned long i = 0; i < vertices.size() / 4; ++i) {
1158 auto aid = vertices[4 * i];
1159 auto bid = vertices[4 * i + 1];
1160 auto cid = vertices[4 * i + 2];
1161 auto did = vertices[4 * i + 3];
1162 mesh.createTetrahedron(mesh.vertex(aid),
1163 mesh.vertex(bid),
1164 mesh.vertex(cid),
1165 mesh.vertex(did));
1166 }
1167}
1168
1170 std::string_view meshName,
1171 std::string_view dataName,
1174{
1175 PRECICE_TRACE(meshName, dataName, vertices.size());
1176 PRECICE_CHECK(_state != State::Finalized, "writeData(...) cannot be called after finalize().");
1177 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.");
1178 PRECICE_REQUIRE_DATA_WRITE(meshName, dataName);
1179 // Inconsistent sizes will be handled below
1180 if (vertices.empty() && values.empty()) {
1181 return;
1182 }
1183
1184 WriteDataContext &context = _accessor->writeDataContext(meshName, dataName);
1185
1186 const auto dataDims = context.getDataDimensions();
1187 const auto expectedDataSize = vertices.size() * dataDims;
1188 PRECICE_CHECK(expectedDataSize == values.size(),
1189 "Input sizes are inconsistent attempting to write {}D data \"{}\" to mesh \"{}\". "
1190 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1191 dataDims, dataName, meshName,
1192 vertices.size(), values.size(), expectedDataSize, dataDims, vertices.size());
1193
1194 // Sizes are correct at this point
1195 PRECICE_VALIDATE_DATA(values.data(), values.size()); // TODO Only take span
1196
1197 if (auto index = context.locateInvalidVertexID(vertices); index) {
1198 PRECICE_ERROR("Cannot write data \"{}\" to mesh \"{}\" due to invalid Vertex ID at vertices[{}]. "
1199 "Please make sure you only use the results from calls to setMeshVertex/Vertices().",
1200 dataName, meshName, *index);
1201 }
1202 Event e{fmt::format("writeData.{}_{}", meshName, dataName), profiling::API};
1203 context.writeValuesIntoDataBuffer(vertices, values);
1204}
1205
1207 std::string_view meshName,
1208 std::string_view dataName,
1210 double relativeReadTime,
1211 ::precice::span<double> values) const
1212{
1213 PRECICE_TRACE(meshName, dataName, vertices.size(), relativeReadTime);
1214 PRECICE_CHECK(_state != State::Constructed, "readData(...) cannot be called before initialize().");
1215 PRECICE_CHECK(_state != State::Finalized, "readData(...) cannot be called after finalize().");
1216 PRECICE_CHECK(math::smallerEquals(relativeReadTime, _couplingScheme->getNextTimeStepMaxSize()), "readData(...) cannot sample data outside of current time window.");
1217 PRECICE_CHECK(relativeReadTime >= 0, "readData(...) cannot sample data before the current time.");
1218 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);
1219
1220 PRECICE_REQUIRE_DATA_READ(meshName, dataName);
1221
1222 PRECICE_CHECK(_meshLock.check(meshName),
1223 "Cannot read from mesh \"{}\" after it has been reset. Please read data before calling resetMesh().",
1224 meshName);
1225
1226 // Inconsistent sizes will be handled below
1227 if (vertices.empty() && values.empty()) {
1228 return;
1229 }
1230
1231 ReadDataContext &context = _accessor->readDataContext(meshName, dataName);
1232 PRECICE_CHECK(context.hasSamples(), "Data \"{}\" cannot be read from mesh \"{}\" as it contains no samples. "
1233 "This is typically a configuration issue of the data flow. "
1234 "Check if the data is correctly exchanged to this participant \"{}\" and mapped to mesh \"{}\".",
1235 dataName, meshName, _accessorName, meshName);
1236
1237 const auto dataDims = context.getDataDimensions();
1238 const auto expectedDataSize = vertices.size() * dataDims;
1239 PRECICE_CHECK(expectedDataSize == values.size(),
1240 "Input/Output sizes are inconsistent attempting to read {}D data \"{}\" from mesh \"{}\". "
1241 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1242 dataDims, dataName, meshName,
1243 vertices.size(), values.size(), expectedDataSize, dataDims, vertices.size());
1244
1245 if (auto index = context.locateInvalidVertexID(vertices); index) {
1246 PRECICE_ERROR("Cannot read data \"{}\" from mesh \"{}\" due to invalid Vertex ID at vertices[{}]. "
1247 "Please make sure you only use the results from calls to setMeshVertex/Vertices().",
1248 dataName, meshName, *index);
1249 }
1250
1251 Event e{fmt::format("readData.{}_{}", meshName, dataName), profiling::API};
1252
1253 double readTime = _couplingScheme->getTime() + relativeReadTime;
1254 context.readValues(vertices, readTime, values);
1255}
1256
1258 std::string_view meshName,
1259 std::string_view dataName,
1261 double relativeReadTime,
1262 ::precice::span<double> values) const
1263{
1265 PRECICE_TRACE(meshName, dataName, coordinates.size(), relativeReadTime);
1266 PRECICE_CHECK(_state != State::Constructed, "mapAndReadData(...) cannot be called before initialize().");
1267 PRECICE_CHECK(_state != State::Finalized, "mapAndReadData(...) cannot be called after finalize().");
1268 PRECICE_CHECK(math::smallerEquals(relativeReadTime, _couplingScheme->getNextTimeStepMaxSize()), "readData(...) cannot sample data outside of current time window.");
1269 PRECICE_CHECK(relativeReadTime >= 0, "mapAndReadData(...) cannot sample data before the current time.");
1270 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);
1271
1272 PRECICE_REQUIRE_DATA_READ(meshName, dataName);
1273 PRECICE_VALIDATE_DATA(coordinates.begin(), coordinates.size());
1274
1275 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1276 "This participant attempteded to map and read data (via \"mapAndReadData\") from mesh \"{0}\", "
1277 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1278 "mapAndReadData({0}, ...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1279 meshName);
1280 // If an access region is required, we have to check its existence
1281 bool requiresBB = requiresUserDefinedAccessRegion(meshName);
1282 PRECICE_CHECK(!requiresBB || (requiresBB && _accessor->receivedMeshContext(meshName).userDefinedAccessRegion),
1283 "The function \"mapAndReadData\" was called on mesh \"{0}\", "
1284 "but no access region was defined although this is necessary for parallel runs. "
1285 "Please define an access region using \"setMeshAccessRegion()\" before calling \"mapAndReadData()\".",
1286 meshName);
1287
1288 PRECICE_CHECK(!_accessor->receivedMeshContext(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. "
1289 "How should the provided data values be read? Please make sure the mesh \"{1}\" is non-empty within the access region.",
1290 dataName, meshName);
1291
1292 ReadDataContext &dataContext = _accessor->readDataContext(meshName, dataName);
1293 PRECICE_CHECK(dataContext.hasJustInTimeMapping(),
1294 "The function \"mapAndReadData\" was called on mesh \"{0}\", but no matching just-in-time mapping was configured. "
1295 "Please define a mapping in read direction from the mesh \"{0}\" and omit the \"to\" attribute from the definition. "
1296 "Example \"<mapping:nearest-neighbor direction=\"read\" from=\"{0}\" constraint=\"consistent\" />",
1297 meshName);
1298
1299 // Inconsistent sizes will be handled below
1300 if (coordinates.empty() && values.empty()) {
1301 return;
1302 }
1303
1304 Event e{fmt::format("mapAndReadData.{}_{}", meshName, dataName), profiling::API};
1305
1306 // Note that meshName refers to a remote mesh
1307 const auto dataDims = dataContext.getDataDimensions();
1308 const auto dim = dataContext.getSpatialDimensions();
1309 const auto nVertices = (coordinates.size() / dim);
1310
1311 // Check that the vertex is actually within the defined access region
1312 _accessor->receivedMeshContext(meshName).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, "writeAndMapData(...) cannot be called before initialize(), because the mesh to map onto hasn't been received yet.");
1335 PRECICE_CHECK(_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.");
1336 PRECICE_REQUIRE_DATA_WRITE(meshName, dataName);
1337
1338 PRECICE_VALIDATE_DATA(coordinates.begin(), coordinates.size());
1339 PRECICE_VALIDATE_DATA(values.data(), values.size());
1340 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1341 "This participant attempteded to map and read data (via \"writeAndMapData\") from mesh \"{0}\", "
1342 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1343 "writeAndMapData({0}, ...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1344 meshName);
1345 // If an access region is required, we have to check its existence
1346 bool requiresBB = requiresUserDefinedAccessRegion(meshName);
1347 PRECICE_CHECK(!requiresBB || (requiresBB && _accessor->receivedMeshContext(meshName).userDefinedAccessRegion),
1348 "The function \"writeAndMapData\" was called on mesh \"{0}\", "
1349 "but no access region was defined although this is necessary for parallel runs. "
1350 "Please define an access region using \"setMeshAccessRegion()\" before calling \"writeAndMapData()\".",
1351 meshName);
1352
1353 WriteDataContext &dataContext = _accessor->writeDataContext(meshName, dataName);
1354 PRECICE_CHECK(dataContext.hasJustInTimeMapping(),
1355 "The function \"writeAndMapData\" was called on mesh \"{0}\", but no matching just-in-time mapping was configured. "
1356 "Please define a mapping in write direction to the mesh \"{0}\" and omit the \"from\" attribute from the definition. "
1357 "Example \"<mapping:nearest-neighbor direction=\"write\" to=\"{0}\" constraint=\"conservative\" />",
1358 meshName);
1359
1360 // Inconsistent sizes will be handled below
1361 if (coordinates.empty() && values.empty()) {
1362 return;
1363 }
1364
1365 Event e{fmt::format("writeAndMapData.{}_{}", meshName, dataName), profiling::API};
1366
1367 // Note that meshName refers here typically to a remote mesh
1368 const auto dataDims = dataContext.getDataDimensions();
1369 const auto dim = dataContext.getSpatialDimensions();
1370 const auto nVertices = (coordinates.size() / dim);
1371 ReceivedMeshContext &context = _accessor->receivedMeshContext(meshName);
1372
1373 // Check that the vertex is actually within the defined access region
1374 _accessor->receivedMeshContext(meshName).checkVerticesInsideAccessRegion(coordinates, dim, "writeAndMapData");
1375
1376 PRECICE_CHECK(nVertices * dataDims == values.size(),
1377 "Input sizes are inconsistent attempting to write {}D data \"{}\" to mesh \"{}\". "
1378 "You passed {} vertex indices and {} data components, but we expected {} data components ({} x {}).",
1379 dataDims, dataName, meshName,
1380 nVertices, values.size(), nVertices * dataDims, dataDims, nVertices);
1381
1382 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. "
1383 "Where should the provided data go? Please make sure the mesh \"{1}\" is non-empty within the access region.",
1384 dataName, meshName);
1385 dataContext.writeAndMapValues(coordinates, values);
1386}
1387
1389 std::string_view meshName,
1390 std::string_view dataName,
1393{
1395
1396 // Asserts and checks
1397 PRECICE_TRACE(meshName, dataName, vertices.size());
1398 PRECICE_CHECK(_state != State::Finalized, "writeGradientData(...) cannot be called after finalize().");
1399 PRECICE_REQUIRE_DATA_WRITE(meshName, dataName);
1400
1401 // Inconsistent sizes will be handled below
1402 if ((vertices.empty() && gradients.empty()) || !requiresGradientDataFor(meshName, dataName)) {
1403 return;
1404 }
1405
1406 // Get the data
1407 WriteDataContext &context = _accessor->writeDataContext(meshName, dataName);
1408
1409 // Check if the Data object of given mesh has been initialized with gradient data
1410 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);
1411
1412 if (auto index = context.locateInvalidVertexID(vertices); index) {
1413 PRECICE_ERROR("Cannot write gradient data \"{}\" to mesh \"{}\" due to invalid Vertex ID at vertices[{}]. "
1414 "Please make sure you only use the results from calls to setMeshVertex/Vertices().",
1415 dataName, meshName, *index);
1416 }
1417
1418 const auto dataDims = context.getDataDimensions();
1419 const auto meshDims = context.getSpatialDimensions();
1420 const auto gradientComponents = meshDims * dataDims;
1421 const auto expectedComponents = vertices.size() * gradientComponents;
1422 PRECICE_CHECK(expectedComponents == gradients.size(),
1423 "Input sizes are inconsistent attempting to write gradient for data \"{}\" to mesh \"{}\". "
1424 "A single gradient/Jacobian for {}D data on a {}D mesh has {} components. "
1425 "You passed {} vertex indices and {} gradient components, but we expected {} gradient components. ",
1426 dataName, meshName,
1427 dataDims, meshDims, gradientComponents,
1428 vertices.size(), gradients.size(), expectedComponents);
1429
1430 PRECICE_VALIDATE_DATA(gradients.data(), gradients.size());
1431
1432 Event e{fmt::format("writeGradientData.{}_{}", meshName, dataName), profiling::API};
1433
1434 context.writeGradientsIntoDataBuffer(vertices, gradients);
1435}
1436
1438 const std::string_view meshName,
1439 ::precice::span<const double> boundingBox) const
1440{
1441 PRECICE_TRACE(meshName, boundingBox.size());
1442 PRECICE_REQUIRE_MESH_USE(meshName);
1443 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1444 "This participant attempteded to set an access region (via \"setMeshAccessRegion\") on mesh \"{0}\", "
1445 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1446 "setMeshAccessRegion(...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1447 meshName);
1448 PRECICE_CHECK(_state != State::Finalized, "setMeshAccessRegion() cannot be called after finalize().");
1449 PRECICE_CHECK(_state != State::Initialized, "setMeshAccessRegion() needs to be called before initialize().");
1450
1451 // Get the related mesh - setMeshAccessRegion only works for received meshes
1452 ReceivedMeshContext &receivedContext = _accessor->receivedMeshContext(meshName);
1453
1454 PRECICE_CHECK(!receivedContext.userDefinedAccessRegion, "A mesh access region was already defined for mesh \"{}\". setMeshAccessRegion may only be called once per mesh.", receivedContext.mesh->getName());
1455 mesh::Mesh &mesh = *receivedContext.mesh;
1456 int dim = mesh.getDimensions();
1457 PRECICE_CHECK(boundingBox.size() == static_cast<unsigned long>(dim) * 2,
1458 "Incorrect amount of bounding box components attempting to set the bounding box of {}D mesh \"{}\" . "
1459 "You passed {} limits, but we expected {} ({}x2).",
1460 dim, meshName, boundingBox.size(), dim * 2, dim);
1461
1462 // Transform bounds into a suitable format
1463 PRECICE_DEBUG("Define bounding box");
1464 std::vector<double> bounds(dim * 2);
1465
1466 for (int d = 0; d < dim; ++d) {
1467 // Check that min is lower or equal to max
1468 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...]");
1469 bounds[2 * d] = boundingBox[2 * d];
1470 bounds[2 * d + 1] = boundingBox[2 * d + 1];
1471 }
1472 // Create a bounding box
1473 receivedContext.userDefinedAccessRegion = std::make_shared<mesh::BoundingBox>(bounds);
1474 // Expand the mesh associated bounding box
1475 mesh.expandBoundingBox(*receivedContext.userDefinedAccessRegion.get());
1476}
1477
1479 const std::string_view meshName,
1481 ::precice::span<double> coordinates) const
1482{
1483 PRECICE_TRACE(meshName, ids.size(), coordinates.size());
1484 PRECICE_REQUIRE_MESH_USE(meshName);
1485 PRECICE_CHECK(_accessor->isMeshReceived(meshName) && _accessor->isDirectAccessAllowed(meshName),
1486 "This participant attempteded to get mesh vertex IDs and coordinates (via \"getMeshVertexIDsAndCoordinates\") from mesh \"{0}\", "
1487 "but mesh \"{0}\" is either not a received mesh or its api access was not enabled in the configuration. "
1488 "getMeshVertexIDsAndCoordinates(...) is only valid for (<receive-mesh name=\"{0}\" ... api-access=\"true\"/>).",
1489 meshName);
1490 // If an access region is required, we have to check its existence
1491 bool requiresBB = requiresUserDefinedAccessRegion(meshName);
1492 PRECICE_CHECK(!requiresBB || (requiresBB && _accessor->receivedMeshContext(meshName).userDefinedAccessRegion),
1493 "The function \"getMeshVertexIDsAndCoordinates\" was called on mesh \"{0}\", "
1494 "but no access region was defined although this is necessary for parallel runs. "
1495 "Please define an access region using \"setMeshAccessRegion()\" before calling \"getMeshVertexIDsAndCoordinates()\".",
1496 meshName);
1497
1498 PRECICE_DEBUG("Get {} mesh vertices with IDs", ids.size());
1499
1500 // Check, if the requested mesh data has already been received. Otherwise, the function call doesn't make any sense
1501 PRECICE_CHECK((_state == State::Initialized) || _accessor->isMeshProvided(meshName),
1502 "initialize() has to be called before accessing data of the received mesh \"{}\" on participant \"{}\".",
1503 meshName, _accessor->getName());
1504
1505 if (ids.empty() && coordinates.empty()) {
1506 return;
1507 }
1508
1509 Event e{fmt::format("getMeshVertexIDsAndCoordinates.{}", meshName), profiling::API};
1510
1511 const MeshContext &context = _accessor->meshContext(meshName);
1512
1513 auto filteredVertices = _accessor->receivedMeshContext(meshName).filterVerticesToLocalAccessRegion(requiresBB);
1514 const auto meshSize = filteredVertices.size();
1515
1516 const mesh::Mesh &mesh = *(context.mesh);
1517 const auto meshDims = mesh.getDimensions();
1518 PRECICE_CHECK(ids.size() == meshSize,
1519 "Output size is incorrect attempting to get vertex ids of {}D mesh \"{}\". "
1520 "You passed {} vertex indices, but we expected {}. "
1521 "Use getMeshVertexSize(\"{}\") to receive the required amount of vertices.",
1522 meshDims, meshName, ids.size(), meshSize, meshName);
1523 const auto expectedCoordinatesSize = static_cast<unsigned long>(meshDims * meshSize);
1524 PRECICE_CHECK(coordinates.size() == expectedCoordinatesSize,
1525 "Output size is incorrect attempting to get vertex coordinates of {}D mesh \"{}\". "
1526 "You passed {} coordinate components, but we expected {} ({}x{}). "
1527 "Use getMeshVertexSize(\"{}\") and getMeshDimensions(\"{}\") to receive the required amount components",
1528 meshDims, meshName, coordinates.size(), expectedCoordinatesSize, meshSize, meshDims, meshName, meshName);
1529
1530 PRECICE_ASSERT(ids.size() <= mesh.nVertices(), "The queried size exceeds the number of available points.");
1531
1532 Eigen::Map<Eigen::MatrixXd> posMatrix{
1533 coordinates.data(), mesh.getDimensions(), static_cast<EIGEN_DEFAULT_DENSE_INDEX_TYPE>(ids.size())};
1534
1535 for (unsigned long i = 0; i < ids.size(); i++) {
1536 auto localID = filteredVertices[i].get().getID();
1537 PRECICE_ASSERT(mesh.isValidVertexID(localID), i, localID);
1538 ids[i] = localID;
1539 posMatrix.col(i) = filteredVertices[i].get().getCoords();
1540 }
1541}
1542
1544{
1545 // sort meshContexts by name, for communication in right order.
1546 std::sort(_accessor->usedMeshContexts().begin(), _accessor->usedMeshContexts().end(),
1547 [](const MeshContextVariant &lhs, const MeshContextVariant &rhs) -> bool {
1548 return getMesh(lhs).getName() < getMesh(rhs).getName();
1549 });
1550
1551 // Provided meshes need their bounding boxes already for the re-partitioning
1552 for (auto &context : _accessor->providedMeshContexts()) {
1553 context.mesh->computeBoundingBox();
1554 }
1555
1556 // Clear mappings for all meshes
1557 for (auto &variant : _accessor->usedMeshContexts()) {
1558 getMeshContext(variant)->clearMappings();
1559 }
1560
1561 // Compare bounding boxes for all meshes
1562 for (const auto &variant : _accessor->usedMeshContexts()) {
1564 }
1565}
1566
1568{
1569 // We need to do this in two loops: First, communicate the mesh and later compute the partition.
1570 // Originally, this was done in one loop. This however gave deadlock if two meshes needed to be communicated cross-wise.
1571 // Both loops need a different sorting
1572
1573 auto &contexts = _accessor->usedMeshContexts();
1574
1575 std::sort(contexts.begin(), contexts.end(),
1576 [](const MeshContextVariant &lhs, const MeshContextVariant &rhs) -> bool {
1577 return getMesh(lhs).getName() < getMesh(rhs).getName();
1578 });
1579
1580 for (const auto &variant : contexts) {
1581 getPartition(variant).communicate();
1582 }
1583
1584 // for two-level initialization, there is also still communication in partition::compute()
1585 // therefore, we cannot resort here.
1586 // @todo this hacky solution should be removed as part of #633
1587 bool resort = true;
1588 for (auto &m2nPair : _m2ns) {
1589 if (m2nPair.second.m2n->usesTwoLevelInitialization()) {
1590 resort = false;
1591 break;
1592 }
1593 }
1594
1595 if (resort) {
1596 // pull provided meshes up front, to have them ready for the decomposition of the received meshes (for the mappings)
1597 std::stable_partition(contexts.begin(), contexts.end(),
1598 [](const MeshContextVariant &variant) -> bool {
1599 return std::holds_alternative<ProvidedMeshContext *>(variant);
1600 });
1601 }
1602
1603 for (const auto &variant : contexts) {
1604 auto &mesh = getMesh(variant);
1605 getPartition(variant).compute();
1606
1607 // Received meshes can only compute their bounding boxes here
1608 if (std::holds_alternative<ReceivedMeshContext *>(variant)) {
1609 mesh.computeBoundingBox();
1610 }
1611
1612 mesh.allocateDataValues();
1613
1614 // Should be relevant for direct mesh access only
1615 const auto requiredSize = mesh.nVertices();
1616 for (auto &context : _accessor->writeDataContexts()) {
1617 if (context.getMeshName() == mesh.getName()) {
1618 context.resizeBufferTo(requiredSize, std::holds_alternative<ReceivedMeshContext *>(variant));
1619 }
1620 }
1621 }
1622}
1623
1624void ParticipantImpl::computeMappings(std::vector<MappingContext> &contexts, const std::string &mappingType)
1625{
1626 PRECICE_TRACE();
1627 bool anyMappingChanged = false;
1628 for (impl::MappingContext &context : contexts) {
1629 if (not context.mapping->hasComputedMapping()) {
1630 PRECICE_INFO_IF(context.configuredWithAliasTag,
1631 "Automatic RBF mapping alias from mesh \"{}\" to mesh \"{}\" in \"{}\" direction resolves to \"{}\" .",
1632 context.mapping->getInputMesh()->getName(), context.mapping->getOutputMesh()->getName(), mappingType, context.mapping->getName());
1633 PRECICE_INFO("Computing \"{}\" mapping from mesh \"{}\" to mesh \"{}\" in \"{}\" direction.",
1634 context.mapping->getName(), context.mapping->getInputMesh()->getName(), context.mapping->getOutputMesh()->getName(), mappingType);
1635 context.mapping->computeMapping();
1636 anyMappingChanged = true;
1637 }
1638 }
1639 if (anyMappingChanged) {
1640 _accessor->initializeMappingDataCache(mappingType);
1641 }
1642}
1643
1645{
1646 PRECICE_TRACE();
1647 if (!_accessor->hasWriteMappings()) {
1648 return;
1649 }
1650
1651 Event e("mapping", profiling::Fundamental);
1652 computeMappings(_accessor->writeMappingContexts(), "write");
1653 for (auto &context : _accessor->writeDataContexts()) {
1654 if (context.hasMapping()) {
1655 PRECICE_DEBUG("Map initial write data \"{}\" from mesh \"{}\"", context.getDataName(), context.getMeshName());
1656 _executedWriteMappings += context.mapData(std::nullopt, true);
1657 }
1658 }
1659}
1660
1661void ParticipantImpl::mapWrittenData(std::optional<double> after)
1662{
1663 PRECICE_TRACE();
1664 if (!_accessor->hasWriteMappings()) {
1665 return;
1666 }
1667
1668 Event e("mapping", profiling::Fundamental);
1669 computeMappings(_accessor->writeMappingContexts(), "write");
1670 for (auto &context : _accessor->writeDataContexts()) {
1671 if (context.hasMapping()) {
1672 PRECICE_DEBUG("Map write data \"{}\" from mesh \"{}\"", context.getDataName(), context.getMeshName());
1673 _executedWriteMappings += context.mapData(after);
1674 }
1675 }
1676}
1677
1678void ParticipantImpl::trimReadMappedData(double startOfTimeWindow, bool isTimeWindowComplete, const cplscheme::ImplicitData &fromData)
1679{
1680 PRECICE_TRACE();
1681 for (auto &context : _accessor->readDataContexts()) {
1682 if (context.hasMapping()) {
1684 // For serial implicit second, we need to discard everything before startOfTimeWindow to preserve the time window start
1685 // For serial implicit first, we need to discard everything as everything is new
1686 // For parallel implicit, we need to discard everything as everything is new
1687 context.clearToDataFor(fromData);
1688 } else {
1689 context.trimToDataAfterFor(fromData, startOfTimeWindow);
1690 }
1691 }
1692 }
1693}
1694
1696{
1697 PRECICE_TRACE();
1698 if (!_accessor->hasReadMappings()) {
1699 return;
1700 }
1701
1702 Event e("mapping", profiling::Fundamental);
1703 computeMappings(_accessor->readMappingContexts(), "read");
1704 for (auto &context : _accessor->readDataContexts()) {
1705 if (context.hasMapping()) {
1706 PRECICE_DEBUG("Map initial read data \"{}\" to mesh \"{}\"", context.getDataName(), context.getMeshName());
1707 // We always ensure that all read data was mapped
1708 _executedReadMappings += context.mapData(std::nullopt, true);
1709 }
1710 }
1711}
1712
1714{
1715 PRECICE_TRACE();
1716 if (!_accessor->hasReadMappings()) {
1717 return;
1718 }
1719
1720 Event e("mapping", profiling::Fundamental);
1721 computeMappings(_accessor->readMappingContexts(), "read");
1722 for (auto &context : _accessor->readDataContexts()) {
1723 if (context.hasMapping()) {
1724 PRECICE_DEBUG("Map read data \"{}\" to mesh \"{}\"", context.getDataName(), context.getMeshName());
1725 // We always ensure that all read data was mapped
1726 _executedReadMappings += context.mapData();
1727 }
1728 }
1729}
1730
1731void ParticipantImpl::performDataActions(const std::set<action::Action::Timing> &timings)
1732{
1733 PRECICE_TRACE();
1734 for (action::PtrAction &action : _accessor->actions()) {
1735 if (timings.find(action->getTiming()) != timings.end()) {
1736 action->performAction();
1737 }
1738 }
1739}
1740
1742{
1743 PRECICE_TRACE();
1744 if (!_accessor->hasExports()) {
1745 return;
1746 }
1747 PRECICE_DEBUG("Handle exports");
1748 profiling::Event e{"handleExports"};
1749
1750 if (timing == ExportTiming::Initial) {
1751 _accessor->exportInitial();
1752 return;
1753 }
1754
1756 exp.timewindow = _couplingScheme->getTimeWindows() - 1;
1758 exp.complete = _couplingScheme->isTimeWindowComplete();
1759 exp.final = !_couplingScheme->isCouplingOngoing();
1760 exp.time = _couplingScheme->getTime();
1761 _accessor->exportIntermediate(exp);
1762}
1763
1765{
1766 PRECICE_TRACE();
1767 for (auto &context : _accessor->writeDataContexts()) {
1768 // reset the buffered data here
1769 context.resetBufferedData();
1770 }
1771}
1772
1775{
1776 const auto &partConfig = *config.getParticipantConfiguration();
1777 PRECICE_CHECK(partConfig.hasParticipant(_accessorName),
1778 "This participant's name, which was specified in the constructor of the preCICE interface as \"{}\", "
1779 "is not defined in the preCICE configuration. "
1780 "Please double-check the correct spelling.",
1782 return partConfig.getParticipant(_accessorName);
1783}
1784
1795
1796void ParticipantImpl::syncTimestep(double computedTimeStepSize)
1797{
1799 Event e("syncTimestep", profiling::Fundamental);
1801 utils::IntraComm::getCommunication()->send(computedTimeStepSize, 0);
1802 } else {
1804 for (Rank secondaryRank : utils::IntraComm::allSecondaryRanks()) {
1805 double dt;
1806 utils::IntraComm::getCommunication()->receive(dt, secondaryRank);
1807 PRECICE_CHECK(math::equals(dt, computedTimeStepSize),
1808 "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 {}.",
1809 secondaryRank, dt, computedTimeStepSize);
1810 }
1811 }
1812}
1813
1815{
1816 PRECICE_DEBUG("Advance coupling scheme");
1817 Event e("advanceCoupling", profiling::Fundamental);
1818 // Orchestrate local and remote mesh changes
1819 std::vector<MeshID> localChanges;
1820
1821 [[maybe_unused]] auto remoteChanges1 = _couplingScheme->firstSynchronization(localChanges);
1822 _couplingScheme->firstExchange();
1823 // Orchestrate remote mesh changes (local ones were handled in the first sync)
1824 [[maybe_unused]] auto remoteChanges2 = _couplingScheme->secondSynchronization();
1825 _couplingScheme->secondExchange();
1826}
1827
1829{
1830 // Optionally apply some final ping-pong to sync solver that run e.g. with a uni-directional coupling
1831 // afterwards close connections
1832 PRECICE_INFO("{} {}communication channels",
1833 (_waitInFinalize ? "Synchronize participants and close" : "Close"),
1834 (close == CloseChannels::Distributed ? "distributed " : ""));
1835 std::string ping = "ping";
1836 std::string pong = "pong";
1837 for (auto &iter : _m2ns) {
1838 auto bm2n = iter.second;
1839 if (!bm2n.m2n->isConnected()) {
1840 PRECICE_DEBUG("Skipping closure of defective connection with {}", bm2n.remoteName);
1841 continue;
1842 }
1844 auto comm = bm2n.m2n->getPrimaryRankCommunication();
1845 PRECICE_DEBUG("Synchronizing primary rank with {}", bm2n.remoteName);
1846 if (bm2n.isRequesting) {
1847 comm->send(ping, 0);
1848 std::string receive = "init";
1849 comm->receive(receive, 0);
1850 PRECICE_ASSERT(receive == pong);
1851 } else {
1852 std::string receive = "init";
1853 comm->receive(receive, 0);
1854 PRECICE_ASSERT(receive == ping);
1855 comm->send(pong, 0);
1856 }
1857 }
1858 if (close == CloseChannels::Distributed) {
1859 PRECICE_DEBUG("Closing distributed communication with {}", bm2n.remoteName);
1860 bm2n.m2n->closeDistributedConnections();
1861 } else {
1862 PRECICE_DEBUG("Closing communication with {}", bm2n.remoteName);
1863 bm2n.m2n->closeConnection();
1864 }
1865 }
1866}
1867
1868bool ParticipantImpl::requiresUserDefinedAccessRegion(std::string_view meshName) const
1869{
1870 return _accessor->isMeshReceived(meshName) && utils::IntraComm::isParallel();
1871}
1872
1873const mesh::Mesh &ParticipantImpl::mesh(const std::string &meshName) const
1874{
1875 PRECICE_TRACE(meshName);
1876 return *_accessor->meshContext(meshName).mesh;
1877}
1878
1886
1887// Reinitialization
1888
1890{
1891 PRECICE_TRACE();
1893 Event e("remesh.exchangeLocalMeshChanges", profiling::Synchronize);
1894
1895 // Gather local changes
1896 std::vector<double> localMeshChanges;
1897 for (const auto &variant : _accessor->usedMeshContexts()) {
1898 localMeshChanges.push_back(_meshLock.check(getMesh(variant).getName()) ? 0.0 : 1.0);
1899 }
1900 PRECICE_DEBUG("Mesh changes of rank: {}", localMeshChanges);
1901
1902 // TODO implement int version of allreduceSum
1903 std::vector<double> totalMeshChanges(localMeshChanges.size(), 0.0);
1904 utils::IntraComm::allreduceSum(localMeshChanges, totalMeshChanges);
1905
1906 // Convert the doubles to int
1907 MeshChanges totalMeshChangesInt(totalMeshChanges.begin(), totalMeshChanges.end());
1908 PRECICE_DEBUG("Mesh changes of participant: {}", totalMeshChangesInt);
1909 return totalMeshChangesInt;
1910}
1911
1913{
1914 // Clear stamples where changes were detected
1915 std::size_t i = 0;
1916 for (auto &variant : _accessor->usedMeshContexts()) {
1917 if (totalMeshChanges[i] > 0.0) {
1918 getMesh(variant).clearDataStamples();
1919 }
1920 ++i;
1921 }
1922}
1923
1924bool ParticipantImpl::reinitHandshake(bool requestReinit) const
1925{
1926 PRECICE_TRACE();
1928 Event e("remesh.exchangeRemoteMeshChanges", profiling::Synchronize);
1929
1931 PRECICE_DEBUG("Remeshing is{} required by this participant.", (requestReinit ? "" : " not"));
1932
1933 bool swarmReinitRequired = requestReinit;
1934 for (auto &iter : _m2ns) {
1935 PRECICE_DEBUG("Coordinating remeshing with {}", iter.first);
1936 bool received = false;
1937 auto &comm = *iter.second.m2n->getPrimaryRankCommunication();
1938 if (iter.second.isRequesting) {
1939 comm.send(requestReinit, 0);
1940 comm.receive(received, 0);
1941 } else {
1942 comm.receive(received, 0);
1943 comm.send(requestReinit, 0);
1944 }
1945 swarmReinitRequired |= received;
1946 }
1947 PRECICE_DEBUG("Coordinated that overall{} remeshing is required.", (swarmReinitRequired ? "" : " no"));
1948
1949 utils::IntraComm::broadcast(swarmReinitRequired);
1950 return swarmReinitRequired;
1951 } else {
1952 bool swarmReinitRequired = false;
1953 utils::IntraComm::broadcast(swarmReinitRequired);
1954 return swarmReinitRequired;
1955 }
1956}
1957
1958void ParticipantImpl::startProfilingSection(std::string_view sectionName)
1959{
1960 PRECICE_CHECK(std::find(sectionName.begin(), sectionName.end(), '/') == sectionName.end(),
1961 "The provided section name \"{}\" may not contain a forward-slash \"/\"",
1962 sectionName);
1963 _userEvents.emplace_back(sectionName, profiling::Fundamental);
1964}
1965
1967{
1968 PRECICE_CHECK(!_userEvents.empty(), "There is no user-started event to stop.");
1969 _userEvents.pop_back();
1970}
1971
1972} // 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
#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)
#define PRECICE_ASSERT(...)
Definition assertion.hpp:85
int getDimensions() const
Definition Mesh.cpp:99
void clear()
Removes all mesh elements and data values (does not remove data or the bounding boxes).
Definition Mesh.cpp:280
const std::string & getName() const
Returns the name of the mesh, as set in the config file.
Definition Mesh.cpp:242
bool empty() const
Does the mesh contain any vertices?
Definition Mesh.hpp:88
Main class for preCICE XML configuration tree.
@ 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 Device.cpp:35
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
std::deque< profiling::Event > _userEvents
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.
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.
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)
std::vector< int > MeshChanges
How many ranks have changed each used 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 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.
Container and creator for meshes.
Definition Mesh.hpp:38
void clearDataStamples()
Clears all data stamples.
Definition Mesh.cpp:304
Vertex of a mesh.
Definition Vertex.hpp:16
virtual void compareBoundingBoxes()=0
Intersections between bounding boxes around each rank are computed.
virtual void compute()=0
The partition is computed, i.e. the mesh re-partitioned if required and all data structures are set u...
virtual void communicate()=0
The mesh is communicated between both primary ranks (if required)
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:63
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
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:59
contains actions to modify exchanged data.
Definition Action.hpp:6
std::unique_ptr< Action > PtrAction
std::shared_ptr< WatchPoint > PtrWatchPoint
std::string errorInvalidVertexID(int vid)
std::variant< ProvidedMeshContext *, ReceivedMeshContext * > MeshContextVariant
Type alias for variant holding either provided or received mesh context pointers.
static constexpr auto errorInvalidVertexIDRange
std::shared_ptr< ParticipantState > PtrParticipant
MeshContext * getMeshContext(const MeshContextVariant &variant)
Helper to extract base MeshContext pointer from variant.
partition::Partition & getPartition(const MeshContextVariant &variant)
Helper to extract partition from variant.
mesh::Mesh & getMesh(const MeshContextVariant &variant)
Helper to extract mesh from variant.
std::shared_ptr< WatchIntegral > PtrWatchIntegral
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
static constexpr Group API
Convenience instance of the Cat::API.
Definition Event.hpp:25
static constexpr SynchronizeTag Synchronize
Convenience instance of the SynchronizeTag.
Definition Event.hpp:28
static constexpr Group Fundamental
Convenience instance of the Cat::Fundamental.
Definition Event.hpp:22
contains the time interpolation logic.
Definition Sample.hpp:8
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.
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
Holds a data mapping and related information.
mesh::PtrMesh mesh
Mesh holding the geometry data structure.
mapping::Mapping::MeshRequirement meshRequirement
Determines which mesh type has to be provided by the accessor.
Context for a mesh provided by this participant.
Context for a mesh received from another participant.
std::shared_ptr< mesh::BoundingBox > userDefinedAccessRegion
Tightly coupled to the parameters of Participant()
Definition XMLTag.hpp:21