UPstreamWrappingTemplates.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2012-2015 OpenFOAM Foundation
9  Copyright (C) 2019-2022 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "UPstreamWrapping.H"
31 #include "PstreamGlobals.H"
32 
33 // * * * * * * * * * * * * * * * Global Functions * * * * * * * * * * * * * //
34 
35 template<class Type>
37 (
38  Type* values,
39  int count,
40  MPI_Datatype datatype,
41  const label comm
42 )
43 {
44  if (!UPstream::parRun())
45  {
46  return;
47  }
48 
49  profilingPstream::beginTiming();
50 
51  // const int retval =
52  MPI_Bcast
53  (
54  values,
55  count,
56  datatype,
57  0, // (root rank) == UPstream::masterNo()
59  );
60 
61  profilingPstream::addBroadcastTime();
62 }
63 
64 
65 template<class Type>
67 (
68  Type* values,
69  int count,
70  MPI_Datatype datatype,
71  MPI_Op optype,
72  const label comm
73 )
74 {
75  if (!UPstream::parRun())
76  {
77  return;
78  }
79 
80  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
81  {
82  Pout<< "** reducing:";
83  if (count == 1)
84  {
85  Pout<< (*values);
86  }
87  else
88  {
89  Pout<< UList<Type>(values, count);
90  }
91  Pout<< " with comm:" << comm
92  << " warnComm:" << UPstream::warnComm << endl;
93  error::printStack(Pout);
94  }
95 
96  profilingPstream::beginTiming();
97 
98  // const int retval =
99  MPI_Reduce
100  (
101  MPI_IN_PLACE, // recv is also send
102  values,
103  count,
104  datatype,
105  optype,
106  0, // (root rank) == UPstream::masterNo()
108  );
110  profilingPstream::addReduceTime();
111 }
112 
113 
114 template<class Type>
116 (
117  Type* values,
118  int count,
119  MPI_Datatype datatype,
120  MPI_Op optype,
121  const label comm,
122  label* requestID
123 )
124 {
125  if (!UPstream::parRun())
126  {
127  return;
128  }
129 
130  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
131  {
132  if (requestID != nullptr)
133  {
134  Pout<< "** MPI_Iallreduce (non-blocking):";
135  }
136  else
137  {
138  Pout<< "** MPI_Allreduce (blocking):";
139  }
140  if (count == 1)
141  {
142  Pout<< (*values);
143  }
144  else
145  {
146  Pout<< UList<Type>(values, count);
147  }
148  Pout<< " with comm:" << comm
149  << " warnComm:" << UPstream::warnComm << endl;
150  error::printStack(Pout);
151  }
152 
153  profilingPstream::beginTiming();
154 
155  bool handled(false);
156 
157 #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
158  if (requestID != nullptr)
159  {
160  handled = true;
161  MPI_Request request;
162  if
163  (
164  MPI_Iallreduce
165  (
166  MPI_IN_PLACE, // recv is also send
167  values,
168  count,
169  datatype,
170  optype,
172  &request
173  )
174  )
175  {
177  << "MPI_Iallreduce failed for "
178  << UList<Type>(values, count)
180  }
181 
182  *requestID = PstreamGlobals::push_request(request);
183  }
184 #endif
185 
186  if (!handled)
187  {
188  if (requestID != nullptr)
189  {
190  *requestID = -1;
191  }
192  if
193  (
194  MPI_Allreduce
195  (
196  MPI_IN_PLACE, // recv is also send
197  values,
198  count,
199  datatype,
200  optype,
202  )
203  )
204  {
206  << "MPI_Allreduce failed for "
207  << UList<Type>(values, count)
209  }
210  }
212  profilingPstream::addReduceTime();
213 }
214 
215 
216 template<class Type>
218 (
219  const UList<Type>& sendData,
220  UList<Type>& recvData,
221  MPI_Datatype datatype,
222  const label comm,
223  label* requestID
224 )
225 {
226  const label np = UPstream::nProcs(comm);
227 
228  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
229  {
230  if (requestID != nullptr)
231  {
232  Pout<< "** MPI_Ialltoall (non-blocking):";
233  }
234  else
235  {
236  Pout<< "** MPI_Alltoall (blocking):";
237  }
238  Pout<< " np:" << np
239  << " sendData:" << sendData.size()
240  << " with comm:" << comm
241  << " warnComm:" << UPstream::warnComm
242  << endl;
243  error::printStack(Pout);
244  }
245 
246  if (sendData.size() != np || recvData.size() != np)
247  {
249  << "Have " << np << " ranks, but size of sendData:"
250  << sendData.size() << " or recvData:" << recvData.size()
251  << " is different!"
253  }
254 
255  if (!UPstream::parRun())
256  {
257  recvData.deepCopy(sendData);
258  return;
259  }
260 
261  profilingPstream::beginTiming();
262 
263  bool handled(false);
264 
265  #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
266  if (requestID != nullptr)
267  {
268  handled = true;
269  MPI_Request request;
270  if
271  (
272  MPI_Ialltoall
273  (
274  // NOTE: const_cast is a temporary hack for
275  // backward-compatibility with versions of OpenMPI < 1.7.4
276  const_cast<Type*>(sendData.cdata()),
277  1, // one element per rank
278  datatype,
279  recvData.data(),
280  1, // one element per rank
281  datatype,
283  &request
284  )
285  )
286  {
288  << "MPI_Ialltoall [comm: " << comm << "] failed."
289  << " For " << sendData
291  }
292 
293  *requestID = PstreamGlobals::push_request(request);
294  }
295 #endif
296 
297 
298  if (!handled)
299  {
300  if (requestID != nullptr)
301  {
302  *requestID = -1;
303  }
304  if
305  (
306  MPI_Alltoall
307  (
308  // NOTE: const_cast is a temporary hack for
309  // backward-compatibility with versions of OpenMPI < 1.7.4
310  const_cast<Type*>(sendData.cdata()),
311  1, // one element per rank
312  datatype,
313  recvData.data(),
314  1, // one element per rank
315  datatype,
317  )
318  )
319  {
321  << "MPI_Alltoall [comm: " << comm << "] failed."
322  << " For " << sendData
324  }
325  }
327  profilingPstream::addAllToAllTime();
328 }
329 
330 
331 template<class Type>
333 (
334  const Type* sendData,
335  const UList<int>& sendCounts,
336  const UList<int>& sendOffsets,
337 
338  Type* recvData,
339  const UList<int>& recvCounts,
340  const UList<int>& recvOffsets,
341 
342  MPI_Datatype datatype,
343  const label comm,
344  label* requestID
345 )
346 {
347  const label np = UPstream::nProcs(comm);
348 
349  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
350  {
351  if (requestID != nullptr)
352  {
353  Pout<< "** MPI_Ialltoallv (non-blocking):";
354  }
355  else
356  {
357  Pout<< "** MPI_Alltoallv (blocking):";
358  }
359  Pout<< " sendCounts:" << sendCounts
360  << " sendOffsets:" << sendOffsets
361  << " with comm:" << comm
362  << " warnComm:" << UPstream::warnComm
363  << endl;
364  error::printStack(Pout);
365  }
366 
367  if
368  (
369  (sendCounts.size() != np || sendOffsets.size() < np)
370  || (recvCounts.size() != np || recvOffsets.size() < np)
371  )
372  {
374  << "Have " << np << " ranks, but sendCounts:" << sendCounts.size()
375  << ", sendOffsets:" << sendOffsets.size()
376  << ", recvCounts:" << recvCounts.size()
377  << " or recvOffsets:" << recvOffsets.size()
378  << " is different!"
380  }
381 
382  if (!UPstream::parRun())
383  {
384  if (recvCounts[0] != sendCounts[0])
385  {
387  << "Bytes to send " << sendCounts[0]
388  << " does not equal bytes to receive " << recvCounts[0]
390  }
391  std::memmove
392  (
393  recvData,
394  (sendData + sendOffsets[0]),
395  recvCounts[0]*sizeof(Type)
396  );
397  return;
398  }
399 
400  profilingPstream::beginTiming();
401 
402  bool handled(false);
403 
404 #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
405  if (requestID != nullptr)
406  {
407  handled = true;
408  MPI_Request request;
409 
410  if
411  (
412  MPI_Ialltoallv
413  (
414  const_cast<Type*>(sendData),
415  const_cast<int*>(sendCounts.cdata()),
416  const_cast<int*>(sendOffsets.cdata()),
417  datatype,
418  recvData,
419  const_cast<int*>(recvCounts.cdata()),
420  const_cast<int*>(recvOffsets.cdata()),
421  datatype,
423  &request
424  )
425  )
426  {
428  << "MPI_Ialltoallv [comm: " << comm << "] failed."
429  << " For sendCounts " << sendCounts
430  << " recvCounts " << recvCounts
432  }
433 
434  *requestID = PstreamGlobals::push_request(request);
435  }
436 #endif
437 
438  if (!handled)
439  {
440  if (requestID != nullptr)
441  {
442  *requestID = -1;
443  }
444  if
445  (
446  MPI_Alltoallv
447  (
448  const_cast<Type*>(sendData),
449  const_cast<int*>(sendCounts.cdata()),
450  const_cast<int*>(sendOffsets.cdata()),
451  datatype,
452  recvData,
453  const_cast<int*>(recvCounts.cdata()),
454  const_cast<int*>(recvOffsets.cdata()),
455  datatype,
457  )
458  )
459  {
461  << "MPI_Alltoallv [comm: " << comm << "] failed."
462  << " For sendCounts " << sendCounts
463  << " recvCounts " << recvCounts
465  }
466  }
468  profilingPstream::addAllToAllTime();
469 }
470 
471 
472 template<class Type>
474 (
475  const Type* sendData,
476  int sendCount,
477 
478  Type* recvData,
479  int recvCount,
480 
481  MPI_Datatype datatype,
482  const label comm,
483  label* requestID
484 )
485 {
486  if (!UPstream::parRun())
487  {
488  std::memmove(recvData, sendData, recvCount*sizeof(Type));
489  return;
490  }
491 
492  const label np = UPstream::nProcs(comm);
493 
494  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
495  {
496  if (requestID != nullptr)
497  {
498  Pout<< "** MPI_Igather (non-blocking):";
499  }
500  else
501  {
502  Pout<< "** MPI_Gather (blocking):";
503  }
504  Pout<< " np:" << np
505  << " recvCount:" << recvCount
506  << " with comm:" << comm
507  << " warnComm:" << UPstream::warnComm
508  << endl;
509  error::printStack(Pout);
510  }
511 
512  profilingPstream::beginTiming();
513 
514  bool handled(false);
515 
516 #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
517  if (requestID != nullptr)
518  {
519  handled = true;
520  MPI_Request request;
521  if
522  (
523  MPI_Igather
524  (
525  const_cast<Type*>(sendData),
526  sendCount,
527  datatype,
528  recvData,
529  recvCount,
530  datatype,
531  0, // (root rank) == UPstream::masterNo()
533  &request
534  )
535  )
536  {
538  << "MPI_Igather [comm: " << comm << "] failed."
539  << " sendCount " << sendCount
540  << " recvCount " << recvCount
542  }
543 
544  *requestID = PstreamGlobals::push_request(request);
545  }
546 #endif
547 
548  if (!handled)
549  {
550  if (requestID != nullptr)
551  {
552  *requestID = -1;
553  }
554  if
555  (
556  MPI_Gather
557  (
558  const_cast<Type*>(sendData),
559  sendCount,
560  datatype,
561  recvData,
562  recvCount,
563  datatype,
564  0, // (root rank) == UPstream::masterNo()
566  )
567  )
568  {
570  << "MPI_Gather [comm: " << comm << "] failed."
571  << " sendCount " << sendCount
572  << " recvCount " << recvCount
574  }
575  }
577  profilingPstream::addGatherTime();
578 }
579 
580 
581 template<class Type>
583 (
584  const Type* sendData,
585  int sendCount,
586 
587  Type* recvData,
588  int recvCount,
589 
590  MPI_Datatype datatype,
591  const label comm,
592  label* requestID
593 )
594 {
595  if (!UPstream::parRun())
596  {
597  std::memmove(recvData, sendData, recvCount*sizeof(Type));
598  return;
599  }
600 
601  const label np = UPstream::nProcs(comm);
602 
603  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
604  {
605  if (requestID != nullptr)
606  {
607  Pout<< "** MPI_Iscatter (non-blocking):";
608  }
609  else
610  {
611  Pout<< "** MPI_Scatter (blocking):";
612  }
613  Pout<< " np:" << np
614  << " recvCount:" << recvCount
615  << " with comm:" << comm
616  << " warnComm:" << UPstream::warnComm
617  << endl;
618  error::printStack(Pout);
619  }
620 
621  profilingPstream::beginTiming();
622 
623  bool handled(false);
624 
625 #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
626  if (requestID != nullptr)
627  {
628  handled = true;
629  MPI_Request request;
630  if
631  (
632  MPI_Iscatter
633  (
634  const_cast<Type*>(sendData),
635  sendCount,
636  datatype,
637  recvData,
638  recvCount,
639  datatype,
640  0, // (root rank) == UPstream::masterNo()
642  &request
643  )
644  )
645  {
647  << "MPI_Iscatter [comm: " << comm << "] failed."
648  << " sendCount " << sendCount
649  << " recvCount " << recvCount
651  }
652 
653  *requestID = PstreamGlobals::push_request(request);
654  }
655 #endif
656 
657  if (!handled)
658  {
659  if (requestID != nullptr)
660  {
661  *requestID = -1;
662  }
663  if
664  (
665  MPI_Scatter
666  (
667  const_cast<Type*>(sendData),
668  sendCount,
669  datatype,
670  recvData,
671  recvCount,
672  datatype,
673  0, // (root rank) == UPstream::masterNo()
675  )
676  )
677  {
679  << "MPI_Iscatter [comm: " << comm << "] failed."
680  << " sendCount " << sendCount
681  << " recvCount " << recvCount
683  }
684  }
686  profilingPstream::addScatterTime();
687 }
688 
689 
690 template<class Type>
692 (
693  const Type* sendData,
694  int sendCount,
695 
696  Type* recvData,
697  const UList<int>& recvCounts,
698  const UList<int>& recvOffsets,
699 
700  MPI_Datatype datatype,
701  const label comm,
702  label* requestID
703 )
704 {
705  if (!UPstream::parRun())
706  {
707  // recvCounts[0] may be invalid - use sendCount instead
708  std::memmove(recvData, sendData, sendCount*sizeof(Type));
709  return;
710  }
711 
712  const label np = UPstream::nProcs(comm);
713 
714  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
715  {
716  if (requestID != nullptr)
717  {
718  Pout<< "** MPI_Igatherv (non-blocking):";
719  }
720  else
721  {
722  Pout<< "** MPI_Gatherv (blocking):";
723  }
724  Pout<< " np:" << np
725  << " recvCounts:" << recvCounts
726  << " recvOffsets:" << recvOffsets
727  << " with comm:" << comm
728  << " warnComm:" << UPstream::warnComm
729  << endl;
730  error::printStack(Pout);
731  }
732 
733  if
734  (
735  UPstream::master(comm)
736  && (recvCounts.size() != np || recvOffsets.size() < np)
737  )
738  {
739  // Note: allow offsets to be e.g. 1 larger than nProc so we
740  // can easily loop over the result
741 
743  << "Have " << np << " ranks, but recvCounts:" << recvCounts.size()
744  << " or recvOffsets:" << recvOffsets.size()
745  << " is too small!"
747  }
748 
749  profilingPstream::beginTiming();
750 
751  // Ensure send/recv consistency on master
752  if (UPstream::master(comm) && !recvCounts[0])
753  {
754  sendCount = 0;
755  }
756 
757  bool handled(false);
758 
759 #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
760  if (requestID != nullptr)
761  {
762  handled = true;
763  MPI_Request request;
764  if
765  (
766  MPI_Igatherv
767  (
768  const_cast<Type*>(sendData),
769  sendCount,
770  datatype,
771  recvData,
772  const_cast<int*>(recvCounts.cdata()),
773  const_cast<int*>(recvOffsets.cdata()),
774  datatype,
775  0, // (root rank) == UPstream::masterNo()
777  &request
778  )
779  )
780  {
782  << "MPI_Igatherv failed [comm: " << comm << ']'
783  << " sendCount " << sendCount
784  << " recvCounts " << recvCounts
786  }
787 
788  *requestID = PstreamGlobals::push_request(request);
789  }
790 #endif
791 
792  if (!handled)
793  {
794  if (requestID != nullptr)
795  {
796  *requestID = -1;
797  }
798  if
799  (
800  MPI_Gatherv
801  (
802  const_cast<Type*>(sendData),
803  sendCount,
804  datatype,
805  recvData,
806  const_cast<int*>(recvCounts.cdata()),
807  const_cast<int*>(recvOffsets.cdata()),
808  datatype,
809  0, // (root rank) == UPstream::masterNo()
811  )
812  )
813  {
815  << "MPI_Gatherv failed [comm: " << comm << ']'
816  << " sendCount " << sendCount
817  << " recvCounts " << recvCounts
819  }
820  }
822  profilingPstream::addGatherTime();
823 }
824 
825 
826 template<class Type>
828 (
829  const Type* sendData,
830  const UList<int>& sendCounts,
831  const UList<int>& sendOffsets,
832 
833  Type* recvData,
834  int recvCount,
835 
836  MPI_Datatype datatype,
837  const label comm,
838  label* requestID
839 )
840 {
841  if (!UPstream::parRun())
842  {
843  std::memmove(recvData, sendData, recvCount*sizeof(Type));
844  return;
845  }
846 
847  const label np = UPstream::nProcs(comm);
848 
849  if (UPstream::warnComm != -1 && comm != UPstream::warnComm)
850  {
851  if (requestID != nullptr)
852  {
853  Pout<< "** MPI_Iscatterv (non-blocking):";
854  }
855  else
856  {
857  Pout<< "** MPI_Scatterv (blocking):";
858  }
859  Pout<< " np:" << np
860  << " sendCounts:" << sendCounts
861  << " sendOffsets:" << sendOffsets
862  << " with comm:" << comm
863  << " warnComm:" << UPstream::warnComm
864  << endl;
865  error::printStack(Pout);
866  }
867 
868  if
869  (
870  UPstream::master(comm)
871  && (sendCounts.size() != np || sendOffsets.size() < np)
872  )
873  {
874  // Note: allow offsets to be e.g. 1 larger than nProc so we
875  // can easily loop over the result
876 
878  << "Have " << np << " ranks, but sendCounts:" << sendCounts.size()
879  << " or sendOffsets:" << sendOffsets.size()
880  << " is too small!"
882  }
883 
884  profilingPstream::beginTiming();
885 
886  bool handled(false);
887 
888 #if defined(MPI_VERSION) && (MPI_VERSION >= 3)
889  if (requestID != nullptr)
890  {
891  handled = true;
892  MPI_Request request;
893  if
894  (
895  MPI_Iscatterv
896  (
897  const_cast<Type*>(sendData),
898  const_cast<int*>(sendCounts.cdata()),
899  const_cast<int*>(sendOffsets.cdata()),
900  datatype,
901  recvData,
902  recvCount,
903  datatype,
904  0, // (root rank) == UPstream::masterNo()
906  &request
907  )
908  )
909  {
911  << "MPI_Iscatterv [comm: " << comm << "] failed."
912  << " sendCounts " << sendCounts
913  << " sendOffsets " << sendOffsets
915  }
916 
917  *requestID = PstreamGlobals::push_request(request);
918  }
919 #endif
920 
921  if (!handled)
922  {
923  if (requestID != nullptr)
924  {
925  *requestID = -1;
926  }
927  if
928  (
929  MPI_Scatterv
930  (
931  const_cast<Type*>(sendData),
932  const_cast<int*>(sendCounts.cdata()),
933  const_cast<int*>(sendOffsets.cdata()),
934  datatype,
935  recvData,
936  recvCount,
937  datatype,
938  0, // (root rank) == UPstream::masterNo()
940  )
941  )
942  {
944  << "MPI_Scatterv [comm: " << comm << "] failed."
945  << " sendCounts " << sendCounts
946  << " sendOffsets " << sendOffsets
948  }
949  }
950 
951  profilingPstream::addScatterTime();
952 }
953 
954 
955 // ************************************************************************* //
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:118
void scatterv(const Type *sendData, const UList< int > &sendCounts, const UList< int > &sendOffsets, Type *recvData, int recvCount, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
error FatalError
Error stream (stdout output on all processes), with additional &#39;FOAM FATAL ERROR&#39; header text and sta...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:578
T * data() noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:230
void gatherv(const Type *sendData, int sendCount, Type *recvData, const UList< int > &recvCounts, const UList< int > &recvOffsets, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
DynamicList< MPI_Comm > MPICommunicators_
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
void allReduce(Type *values, int count, MPI_Datatype datatype, MPI_Op optype, const label comm, label *requestID=nullptr)
void broadcast0(Type *values, int count, MPI_Datatype datatype, const label comm)
List< T > values(const HashTable< T, Key, Hash > &tbl, const bool doSort=false)
List of values from HashTable, optionally sorted.
Definition: HashOps.H:164
unsigned int count(const UList< bool > &bools, const bool val=true)
Count number of &#39;true&#39; entries.
Definition: BitOps.H:73
void gather(const Type *sendData, int sendCount, Type *recvData, int recvCount, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
errorManip< error > abort(error &err)
Definition: errorManip.H:139
void allToAll(const UList< Type > &sendData, UList< Type > &recvData, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
label push_request(MPI_Request request)
Reuse previously freed request locations or push request onto list of outstanding requests...
void allToAllv(const Type *sendData, const UList< int > &sendCounts, const UList< int > &sendOffsets, Type *recvData, const UList< int > &recvCounts, const UList< int > &recvOffsets, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
const T * cdata() const noexcept
Return pointer to the underlying array serving as data storage.
Definition: UListI.H:223
void deepCopy(const UList< T > &list)
Copy elements of the given UList. Sizes must match!
Definition: UList.C:99
void reduce0(Type *values, int count, MPI_Datatype datatype, MPI_Op optype, const label comm)
void scatter(const Type *sendData, int sendCount, Type *recvData, int recvCount, MPI_Datatype datatype, const label comm, label *requestID=nullptr)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.