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