Simplify mdrun tMPI thread spawning
[gromacs.git] / src / external / thread_mpi / src / tmpi_init.cpp
blobe2576ef0bb3dca9ad4aad370e7bf5aac22eab8ef
1 /*
2 This source code file is part of thread_mpi.
3 Written by Sander Pronk, Erik Lindahl, and possibly others.
5 Copyright (c) 2009,2018, Sander Pronk, Erik Lindahl.
6 All rights reserved.
8 Redistribution and use in source and binary forms, with or without
9 modification, are permitted provided that the following conditions are met:
10 1) Redistributions of source code must retain the above copyright
11 notice, this list of conditions and the following disclaimer.
12 2) Redistributions in binary form must reproduce the above copyright
13 notice, this list of conditions and the following disclaimer in the
14 documentation and/or other materials provided with the distribution.
15 3) Neither the name of the copyright holders nor the
16 names of its contributors may be used to endorse or promote products
17 derived from this software without specific prior written permission.
19 THIS SOFTWARE IS PROVIDED BY US ''AS IS'' AND ANY
20 EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
21 WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
22 DISCLAIMED. IN NO EVENT SHALL WE BE LIABLE FOR ANY
23 DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
24 (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
25 LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND
26 ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
27 (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
28 SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
30 If you want to redistribute modifications, please consider that
31 scientific software is very special. Version control is crucial -
32 bugs must be traceable. We will be happy to consider code for
33 inclusion in the official distribution, but derived work should not
34 be called official thread_mpi. Details are found in the README & COPYING
35 files.
39 #ifdef HAVE_TMPI_CONFIG_H
40 #include "tmpi_config.h"
41 #endif
43 #ifdef HAVE_CONFIG_H
44 #include "config.h"
45 #endif
47 #ifdef HAVE_UNISTD_H
48 #include <unistd.h>
49 #endif
51 #include <errno.h>
52 #include <stdlib.h>
53 #include <stdio.h>
54 #include <string.h>
55 #if !(defined( _WIN32 ) || defined( _WIN64 ) )
56 #include <sys/time.h>
58 #endif
60 #include <cassert>
62 #include "impl.h"
64 #ifdef TMPI_TRACE
65 #include <stdarg.h>
66 #endif
73 /* there are a few global variables that maintain information about the
74 running threads. Some are defined by the MPI standard: */
75 /* TMPI_COMM_WORLD is in tmpi_malloc.c due to technical reasons */
76 tMPI_Group TMPI_GROUP_EMPTY = NULL;
79 /* the threads themselves (tmpi_comm only contains lists of pointers to this
80 structure */
81 struct tmpi_thread *threads = NULL;
82 int Nthreads = 0;
84 /* thread info */
85 tMPI_Thread_key_t id_key; /* the key to get the thread id */
89 /* whether MPI has finalized (we need this to distinguish pre-inited from
90 post-finalized states */
91 static tmpi_bool tmpi_finalized = FALSE;
93 /* misc. global information about MPI */
94 struct tmpi_global *tmpi_global = NULL;
101 /* start N threads with argc, argv (used by tMPI_Init)*/
102 static
103 int tMPI_Start_threads(tmpi_bool main_returns, int N,
104 tMPI_Affinity_strategy aff_strategy,
105 int *argc, char ***argv,
106 void (*start_fn)(const void*), const void *start_arg,
107 int (*start_fn_main)(int, char**));
109 /* starter function for threads; takes a void pointer to a
110 struct tmpi_starter_, which calls main() if tmpi_start_.fn == NULL */
111 static void* tMPI_Thread_starter(void *arg);
113 /* allocate and initialize the data associated with a thread structure */
114 static int tMPI_Thread_init(struct tmpi_thread *th);
115 /* deallocate the data associated with a thread structure */
116 static void tMPI_Thread_destroy(struct tmpi_thread *th);
121 #ifdef TMPI_TRACE
122 void tMPI_Trace_print(const char *fmt, ...)
124 va_list argp;
125 struct tmpi_thread * th = NULL;
126 static tMPI_Thread_mutex_t mtx = TMPI_THREAD_MUTEX_INITIALIZER;
128 /* don't check for errors during trace */
129 tMPI_Thread_mutex_lock(&mtx);
130 if (threads)
132 th = tMPI_Get_current();
133 printf("THREAD %02d: ", (int)(th-threads));
135 else
137 printf("THREAD main: ");
139 va_start(argp, fmt);
140 vprintf(fmt, argp);
141 printf("\n");
142 fflush(stdout);
143 va_end(argp);
144 tMPI_Thread_mutex_unlock(&mtx);
146 #endif
149 tmpi_bool tMPI_Is_master(void)
151 /* if there are no other threads, we're the main thread */
152 if ( (!TMPI_COMM_WORLD) || TMPI_COMM_WORLD->grp.N == 0)
154 return TRUE;
157 /* otherwise we know this through thread specific data: */
158 /* whether the thread pointer points to the head of the threads array */
159 return (tmpi_bool)(tMPI_Get_current() == threads);
162 tMPI_Comm tMPI_Get_comm_self(void)
164 struct tmpi_thread* th = tMPI_Get_current();
165 return th->self_comm;
169 int tMPI_Get_N(int *argc, char ***argv, const char *optname, int *nthreads)
171 int i;
172 int ret = TMPI_SUCCESS;
174 *nthreads = 0;
175 if (!optname)
177 i = 0;
179 else
181 for (i = 1; i < *argc; i++)
183 if (strcmp(optname, (*argv)[i]) == 0)
185 break;
189 if (i+1 < (*argc))
191 /* the number of processes is an argument */
192 char *end;
193 *nthreads = strtol((*argv)[i+1], &end, 10);
194 if (!end || (*end != 0) )
196 *nthreads = 0;
197 ret = TMPI_FAILURE;
200 if (*nthreads < 1)
202 int nth = tMPI_Thread_get_hw_number();
204 if (nth < 1)
206 nth = 1; /* make sure it's at least 1 */
208 *nthreads = nth;
211 return ret;
214 static int tMPI_Thread_init(struct tmpi_thread *th)
216 int ret;
217 int N_envelopes = (Nthreads+1)*N_EV_ALLOC;
218 int N_send_envelopes = N_EV_ALLOC;
219 int N_reqs = (Nthreads+1)*N_EV_ALLOC;
220 int i;
222 /* we set our thread id, as a thread-specific piece of global data. */
223 ret = tMPI_Thread_setspecific(id_key, th);
224 if (ret != 0)
226 return ret;
229 /* allocate comm.self */
230 ret = tMPI_Comm_alloc( &(th->self_comm), TMPI_COMM_WORLD, 1);
231 if (ret != TMPI_SUCCESS)
233 return ret;
235 th->self_comm->grp.peers[0] = th;
237 /* allocate envelopes */
238 ret = tMPI_Free_env_list_init( &(th->envelopes), N_envelopes );
239 if (ret != TMPI_SUCCESS)
241 return ret;
243 /* recv list */
244 ret = tMPI_Recv_env_list_init( &(th->evr));
245 if (ret != TMPI_SUCCESS)
247 return ret;
249 /* send lists */
250 th->evs = (struct send_envelope_list*)tMPI_Malloc(
251 sizeof(struct send_envelope_list)*Nthreads);
252 if (th->evs == NULL)
254 return TMPI_ERR_NO_MEM;
256 for (i = 0; i < Nthreads; i++)
258 ret = tMPI_Send_env_list_init( &(th->evs[i]), N_send_envelopes);
259 if (ret != TMPI_SUCCESS)
261 return ret;
265 tMPI_Atomic_set( &(th->ev_outgoing_received), 0);
267 tMPI_Event_init( &(th->p2p_event) );
269 /* allocate requests */
270 ret = tMPI_Req_list_init(&(th->rql), N_reqs);
271 if (ret != TMPI_SUCCESS)
273 return ret;
277 #ifdef USE_COLLECTIVE_COPY_BUFFER
278 /* allcate copy_buffer list */
279 ret = tMPI_Copy_buffer_list_init(&(th->cbl_multi),
280 (Nthreads+1)*(N_COLL_ENV+1),
281 Nthreads*COPY_BUFFER_SIZE);
282 if (ret != TMPI_SUCCESS)
284 return ret;
286 #endif
288 #ifdef TMPI_PROFILE
289 ret = tMPI_Profile_init(&(th->profile));
290 if (ret != TMPI_SUCCESS)
292 return ret;
294 #endif
295 /* now wait for all other threads to come on line, before we
296 start the MPI program */
297 ret = tMPI_Thread_barrier_wait( &(tmpi_global->barrier) );
298 if (ret != 0)
300 return ret;;
302 return ret;
306 static void tMPI_Thread_destroy(struct tmpi_thread *th)
308 int i;
310 tMPI_Recv_env_list_destroy( &(th->evr));
311 for (i = 0; i < Nthreads; i++)
313 tMPI_Send_env_list_destroy( &(th->evs[i]));
315 free(th->evs);
316 tMPI_Free_env_list_destroy( &(th->envelopes) );
317 tMPI_Event_destroy( &(th->p2p_event) );
318 tMPI_Req_list_destroy( &(th->rql) );
320 #ifdef USE_COLLECTIVE_COPY_BUFFER
321 tMPI_Copy_buffer_list_destroy(&(th->cbl_multi));
322 #endif
324 for (i = 0; i < th->argc; i++)
326 free(th->argv[i]);
330 static int tMPI_Global_init(struct tmpi_global *g, int Nthreads)
332 int ret;
334 g->usertypes = NULL;
335 g->N_usertypes = 0;
336 g->Nalloc_usertypes = 0;
337 ret = tMPI_Thread_mutex_init(&(g->timer_mutex));
338 if (ret != 0)
340 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
342 tMPI_Spinlock_init(&(g->datatype_lock));
344 ret = tMPI_Thread_barrier_init( &(g->barrier), Nthreads);
345 if (ret != 0)
347 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
350 ret = tMPI_Thread_mutex_init(&(g->comm_link_lock));
351 if (ret != 0)
353 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
357 #if !(defined( _WIN32 ) || defined( _WIN64 ) )
358 /* the time at initialization. */
359 gettimeofday( &(g->timer_init), NULL);
360 #else
361 /* the time at initialization. */
362 g->timer_init = GetTickCount();
363 #endif
364 return TMPI_SUCCESS;
367 static void tMPI_Global_destroy(struct tmpi_global *g)
369 tMPI_Thread_barrier_destroy(&(g->barrier));
370 tMPI_Thread_mutex_destroy(&(g->timer_mutex));
371 tMPI_Thread_mutex_destroy(&(g->comm_link_lock));
377 static void* tMPI_Thread_starter(void *arg)
379 int ret;
380 struct tmpi_thread *th = (struct tmpi_thread*)arg;
382 #ifdef TMPI_TRACE
383 tMPI_Trace_print("Created thread nr. %d", (int)(th-threads));
384 #endif
386 ret = tMPI_Thread_init(th);
387 if (ret != TMPI_SUCCESS)
389 return NULL;
392 /* start_fn, start_arg, argc and argv were set by the calling function */
393 if (!th->start_fn)
395 th->start_fn_main(th->argc, th->argv);
397 else
399 th->start_fn(th->start_arg);
400 if (!tmpi_finalized)
402 tMPI_Finalize();
406 return NULL;
410 int tMPI_Start_threads(tmpi_bool main_returns, int N,
411 tMPI_Affinity_strategy aff_strategy,
412 int *argc, char ***argv,
413 void (*start_fn)(const void*), const void *start_arg,
414 int (*start_fn_main)(int, char**))
416 int ret;
417 #ifdef TMPI_TRACE
418 tMPI_Trace_print("tMPI_Start_threads(%d, %d, %d, %d, %d, %p, %p, %p, %p)",
419 main_returns, N, aff_strategy, argc, argv, start_fn,
420 start_arg);
421 #endif
422 if (N > 0)
424 int i;
425 int set_affinity = FALSE;
427 tmpi_finalized = FALSE;
428 Nthreads = N;
430 /* allocate global data */
431 tmpi_global = (struct tmpi_global*)
432 tMPI_Malloc(sizeof(struct tmpi_global));
433 if (tmpi_global == 0)
435 return TMPI_ERR_NO_MEM;
437 ret = tMPI_Global_init(tmpi_global, N);
438 if (ret != TMPI_SUCCESS)
440 return ret;
443 /* allocate world and thread data */
444 threads = (struct tmpi_thread*)
445 tMPI_Malloc(sizeof(struct tmpi_thread)*N);
446 if (threads == NULL)
448 return TMPI_ERR_NO_MEM;
450 ret = tMPI_Comm_alloc(&TMPI_COMM_WORLD, NULL, N);
451 if (ret != TMPI_SUCCESS)
453 return ret;
455 assert(TMPI_COMM_WORLD != nullptr);
456 TMPI_GROUP_EMPTY = tMPI_Group_alloc();
458 if (tMPI_Thread_key_create(&id_key, NULL))
460 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_INIT);
462 for (i = 0; i < N; i++)
464 TMPI_COMM_WORLD->grp.peers[i] = &(threads[i]);
466 /* copy argc, argv */
467 if (argc && argv)
469 int j;
470 threads[i].argc = *argc;
471 threads[i].argv = (char**)tMPI_Malloc(threads[i].argc*
472 sizeof(char*));
473 for (j = 0; j < threads[i].argc; j++)
475 #if !(defined( _WIN32 ) || defined( _WIN64 ) )
476 threads[i].argv[j] = strdup( (*argv)[j] );
477 #else
478 threads[i].argv[j] = _strdup( (*argv)[j] );
479 #endif
482 else
484 threads[i].argc = 0;
485 threads[i].argv = NULL;
487 threads[i].start_fn = start_fn;
488 threads[i].start_fn_main = start_fn_main;
489 threads[i].start_arg = start_arg;
492 /* now check whether to set affinity */
493 if (aff_strategy == TMPI_AFFINITY_ALL_CORES)
495 int nhw = tMPI_Thread_get_hw_number();
496 if ((nhw > 1) && (nhw == N))
498 set_affinity = TRUE;
502 /* set thread 0's properties */
503 threads[0].thread_id = tMPI_Thread_self();
504 if (set_affinity)
506 /* set the main thread's affinity */
507 tMPI_Thread_setaffinity_single(threads[0].thread_id, 0);
510 for (i = 1; i < N; i++) /* zero is the main thread */
512 ret = tMPI_Thread_create(&(threads[i].thread_id),
513 tMPI_Thread_starter,
514 (void*)&(threads[i]) );
516 if (set_affinity)
518 tMPI_Thread_setaffinity_single(threads[i].thread_id, i);
520 if (ret != TMPI_SUCCESS)
522 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_INIT);
525 /* the main thread also runs start_fn if we don't want
526 it to return */
527 if (!main_returns)
529 tMPI_Thread_starter((void*)&(threads[0]));
532 else
534 ret = tMPI_Thread_init(&(threads[0]));
535 if (ret != 0)
537 return ret;
541 return TMPI_SUCCESS;
545 int tMPI_Init(int *argc, char ***argv,
546 int (*start_function)(int, char**))
548 int ret;
549 #ifdef TMPI_TRACE
550 tMPI_Trace_print("tMPI_Init(%p, %p, %p)", argc, argv, start_function);
551 #endif
553 if (TMPI_COMM_WORLD == 0) /* we're the main process */
555 int N = 0;
556 tMPI_Get_N(argc, argv, "-nt", &N);
557 ret = tMPI_Start_threads(TRUE, N, TMPI_AFFINITY_ALL_CORES, argc, argv,
558 NULL, NULL, start_function) != 0;
559 if (ret != 0)
561 return ret;
564 else
566 /* if we're a sub-thread we need don't need to do anyhing, because
567 everything has already been set up by either the main thread,
568 or the thread runner function.*/
570 return TMPI_SUCCESS;
576 int tMPI_Init_fn(int main_thread_returns, int N,
577 tMPI_Affinity_strategy aff_strategy,
578 void (*start_function)(const void*), const void *arg)
580 int ret;
581 #ifdef TMPI_TRACE
582 tMPI_Trace_print("tMPI_Init_fn(%d, %p, %p)", N, start_function, arg);
583 #endif
585 if (N < 1)
587 N = tMPI_Thread_get_hw_number();
588 if (N < 1)
590 N = 1; /*because that's what the fn returns if it doesn't know*/
594 if (TMPI_COMM_WORLD == 0 && N >= 1) /* we're the main process */
596 ret = tMPI_Start_threads(main_thread_returns, N, aff_strategy,
597 0, 0, start_function, arg, NULL);
598 if (ret != 0)
600 return ret;
603 return TMPI_SUCCESS;
606 int tMPI_Initialized(int *flag)
608 #ifdef TMPI_TRACE
609 tMPI_Trace_print("tMPI_Initialized(%p)", flag);
610 #endif
612 *flag = (TMPI_COMM_WORLD && !tmpi_finalized);
614 return TMPI_SUCCESS;
617 int tMPI_Finalize(void)
619 int i;
620 int ret;
621 #ifdef TMPI_TRACE
622 tMPI_Trace_print("tMPI_Finalize()");
623 #endif
624 #ifdef TMPI_DEBUG
625 printf("%5d: tMPI_Finalize called\n", tMPI_This_threadnr());
626 fflush(stdout);
627 #endif
629 #ifdef TMPI_PROFILE
631 struct tmpi_thread *cur = tMPI_Get_current();
633 tMPI_Profile_stop( &(cur->profile) );
634 ret = tMPI_Thread_barrier_wait( &(tmpi_global->barrier) );
635 if (ret != 0)
637 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
640 if (tMPI_Is_master())
642 tMPI_Profiles_summarize(Nthreads, threads);
645 #endif
646 ret = tMPI_Thread_barrier_wait( &(tmpi_global->barrier) );
647 if (ret != 0)
649 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
654 if (tMPI_Is_master())
657 /* we just wait for all threads to finish; the order isn't very
658 relevant, as all threads should arrive at their endpoints soon. */
659 for (i = 1; i < Nthreads; i++)
661 if (tMPI_Thread_join(threads[i].thread_id, NULL))
663 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_FINALIZE);
665 tMPI_Thread_destroy(&(threads[i]));
667 /* at this point, we are the only thread left, so we can
668 destroy the global structures with impunity. */
669 tMPI_Thread_destroy(&(threads[0]));
670 free(threads);
672 tMPI_Thread_key_delete(id_key);
673 /* de-allocate all the comm stuctures. */
675 tMPI_Comm cur;
677 ret = tMPI_Thread_mutex_lock(&(tmpi_global->comm_link_lock));
678 if (ret != 0)
680 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
682 cur = TMPI_COMM_WORLD->next;
683 while (cur && (cur != TMPI_COMM_WORLD) )
685 tMPI_Comm next = cur->next;
686 ret = tMPI_Comm_destroy(cur, FALSE);
687 if (ret != 0)
689 tMPI_Thread_mutex_unlock(&(tmpi_global->comm_link_lock));
690 return ret;
692 cur = next;
694 ret = tMPI_Comm_destroy(TMPI_COMM_WORLD, FALSE);
695 if (ret != 0)
697 tMPI_Thread_mutex_unlock(&(tmpi_global->comm_link_lock));
698 return ret;
700 ret = tMPI_Thread_mutex_unlock(&(tmpi_global->comm_link_lock));
701 if (ret != 0)
703 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_IO);
708 tMPI_Group_free(&TMPI_GROUP_EMPTY);
709 threads = 0;
710 TMPI_COMM_WORLD = NULL;
711 TMPI_GROUP_EMPTY = NULL;
712 Nthreads = 0;
714 /* deallocate the 'global' structure */
715 tMPI_Global_destroy(tmpi_global);
716 free(tmpi_global);
718 tmpi_finalized = TRUE;
720 else
722 tMPI_Thread_exit(0);
724 return TMPI_SUCCESS;
728 int tMPI_Finalized(int *flag)
730 #ifdef TMPI_TRACE
731 tMPI_Trace_print("tMPI_Finalized(%p)", flag);
732 #endif
733 *flag = tmpi_finalized;
735 return TMPI_SUCCESS;
740 int tMPI_Abort(tMPI_Comm comm, int errorcode)
742 #ifdef TMPI_TRACE
743 tMPI_Trace_print("tMPI_Abort(%p, %d)", comm, errorcode);
744 #endif
745 #if 0
746 /* we abort(). This way we can run a debugger on it */
747 fprintf(stderr, "tMPI_Abort called with error code %d", errorcode);
748 if (comm == TMPI_COMM_WORLD)
750 fprintf(stderr, " on TMPI_COMM_WORLD");
752 fprintf(stderr, "\n");
753 fflush(stdout);
755 abort();
756 #else
757 /* we just kill all threads, but not the main process */
759 if (tMPI_Is_master())
761 if (comm == TMPI_COMM_WORLD)
763 fprintf(stderr,
764 "tMPI_Abort called on TMPI_COMM_WORLD main with errorcode=%d\n",
765 errorcode);
767 else
769 fprintf(stderr,
770 "tMPI_Abort called on main thread with errorcode=%d\n",
771 errorcode);
773 fflush(stderr);
774 exit(errorcode);
776 else
778 int *ret;
779 /* kill myself */
780 fprintf(stderr, "tMPI_Abort called with error code %d on thread %d\n",
781 errorcode, tMPI_This_threadnr());
782 fflush(stderr);
783 ret = (int*)malloc(sizeof(int));
784 tMPI_Thread_exit(ret);
786 #endif
787 return TMPI_SUCCESS;
791 int tMPI_Get_processor_name(char *name, int *resultlen)
793 int nr = tMPI_Threadnr(tMPI_Get_current());
794 unsigned int digits = 0;
795 const unsigned int base = 10;
797 #ifdef TMPI_TRACE
798 tMPI_Trace_print("tMPI_Get_processor_name(%p, %p)", name, resultlen);
799 #endif
800 /* we don't want to call sprintf here (it turns out to be not entirely
801 thread-safe on Mac OS X, for example), so we do it our own way: */
803 /* first determine number of digits */
805 int rest = nr;
806 while (rest > 0)
808 rest /= base;
809 digits++;
811 if (digits == 0)
813 digits = 1;
816 #ifndef _MSC_VER
817 strcpy(name, "thread #");
818 #else
819 strncpy_s(name, TMPI_MAX_PROCESSOR_NAME, "thread #", TMPI_MAX_PROCESSOR_NAME);
820 #endif
821 /* now construct the number */
823 size_t len = strlen(name);
824 unsigned int i;
825 int rest = nr;
827 for (i = 0; i < digits; i++)
829 size_t pos = len + (digits-i-1);
830 if (pos < (TMPI_MAX_PROCESSOR_NAME -1) )
832 name[ pos ] = (char)('0' + rest%base);
834 rest /= base;
836 if ( (digits+len) < TMPI_MAX_PROCESSOR_NAME)
838 name[digits + len] = '\0';
840 else
842 name[TMPI_MAX_PROCESSOR_NAME] = '\0';
846 if (resultlen)
848 *resultlen = (int)strlen(name); /* For some reason the MPI standard
849 uses ints instead of size_ts for
850 sizes. */
852 return TMPI_SUCCESS;
859 /* TODO: there must be better ways to do this */
860 double tMPI_Wtime(void)
862 double ret = 0;
864 #ifdef TMPI_TRACE
865 tMPI_Trace_print("tMPI_Wtime()");
866 #endif
868 #if !(defined( _WIN32 ) || defined( _WIN64 ) )
870 struct timeval tv;
871 long int secdiff;
872 int usecdiff;
874 gettimeofday(&tv, NULL);
875 secdiff = tv.tv_sec - tmpi_global->timer_init.tv_sec;
876 usecdiff = tv.tv_usec - tmpi_global->timer_init.tv_usec;
878 ret = (double)secdiff + 1e-6*usecdiff;
880 #else
882 DWORD tv = GetTickCount();
884 /* the windows absolute time GetTickCount() wraps around in ~49 days,
885 so it's safer to always use differences, and assume that our
886 program doesn't run that long.. */
887 ret = 1e-3*((unsigned int)(tv - tmpi_global->timer_init));
889 #endif
890 return ret;
893 double tMPI_Wtick(void)
895 #if !(defined( _WIN32 ) || defined( _WIN64 ) )
896 /* In Unix, we don't really know. Any modern OS should be at least
897 this precise, though */
898 return 1./100.;
899 #else
900 /* According to the Windows documentation, this is about right: */
901 return 1./100.;
902 #endif
905 int tMPI_Get_count(tMPI_Status *status, tMPI_Datatype datatype, int *count)
907 #ifdef TMPI_TRACE
908 tMPI_Trace_print("tMPI_Get_count(%p, %p, %p)", status, datatype, count);
909 #endif
910 if (!status)
912 return tMPI_Error(TMPI_COMM_WORLD, TMPI_ERR_STATUS);
914 *count = (int)(status->transferred/datatype->size);
915 return TMPI_SUCCESS;