#include #include #include #include #include #include #include #include #include #include #include #include #include #include #include "adler32.hh" namespace Carpet { using namespace std; static void CallScheduledFunction (char const * restrict time_and_mode, void * function, cFunctionData * attribute, void * data, Timers::Timer & user_timer); static void SyncGroupsInScheduleBlock (cFunctionData * attribute, cGH * cctkGH, vector const & sync_groups, Timers::Timer & sync_timer); /// Traverse one function on all components of one refinement level /// of one multigrid level. int CallFunction (void * function, ///< the function to call cFunctionData * attribute, ///< attributes of the function void * data) ///< private data for CCTK_CallFunction { DECLARE_CCTK_PARAMETERS; static Timers::Timer total_timer ("CallFunction"); static Timers::Timer user_timer ("thorns"); static Timers::Timer sync_timer ("syncs"); total_timer.start(); cGH * cctkGH = static_cast (data); assert (int (not not attribute->meta) + int (not not attribute->meta_early) + int (not not attribute->meta_late) + int (not not attribute->global) + int (not not attribute->global_early) + int (not not attribute->global_late) + int (not not attribute->level) + int (not not attribute->singlemap) + int (not not attribute->local) <= 1); assert (not not attribute->loop_global + not not attribute->loop_level + not not attribute->loop_singlemap + not not attribute->loop_local <= 1); // Create list of all groups that need to be synchronised vector sync_groups; sync_groups.reserve (attribute->n_SyncGroups); for (int g = 0; g < attribute->n_SyncGroups; g++) { const int group = attribute->SyncGroups[g]; if (CCTK_NumVarsInGroupI (group) > 0) { // don't add empty groups from the list sync_groups.push_back (group); } } if (attribute->meta or attribute->meta_early or attribute->meta_late or is_meta_mode()) { // Convtest operation if ((attribute->meta and do_meta_mode) or (attribute->meta_early and do_early_meta_mode) or (attribute->meta_late and do_late_meta_mode) or is_meta_mode()) { if (attribute->loop_local) { BEGIN_META_MODE(cctkGH) { BEGIN_MGLEVEL_LOOP(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { BEGIN_LOCAL_MAP_LOOP(cctkGH, CCTK_GF) { BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Meta time local mode", function, attribute, data, user_timer); } END_LOCAL_COMPONENT_LOOP; } END_LOCAL_MAP_LOOP; if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } END_REFLEVEL_LOOP; } END_MGLEVEL_LOOP; } END_META_MODE; } else if (attribute->loop_singlemap) { BEGIN_META_MODE(cctkGH) { BEGIN_MGLEVEL_LOOP(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { BEGIN_MAP_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Meta time singlemap mode", function, attribute, data, user_timer); } END_MAP_LOOP; if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } END_REFLEVEL_LOOP; } END_MGLEVEL_LOOP; } END_META_MODE; } else if (attribute->loop_level) { BEGIN_META_MODE(cctkGH) { BEGIN_MGLEVEL_LOOP(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { CallScheduledFunction ("Meta time level mode", function, attribute, data, user_timer); if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } END_REFLEVEL_LOOP; } END_MGLEVEL_LOOP; } END_META_MODE; } else if (attribute->loop_global) { BEGIN_META_MODE(cctkGH) { BEGIN_MGLEVEL_LOOP(cctkGH) { CallScheduledFunction ("Meta time global mode", function, attribute, data, user_timer); if (not sync_groups.empty()) { BEGIN_REFLEVEL_LOOP(cctkGH) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } END_REFLEVEL_LOOP; } } END_MGLEVEL_LOOP; } END_META_MODE; } else { BEGIN_META_MODE(cctkGH) { CallScheduledFunction ("Meta mode", function, attribute, data, user_timer); if (not sync_groups.empty()) { BEGIN_MGLEVEL_LOOP(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } END_REFLEVEL_LOOP; } END_MGLEVEL_LOOP; } } END_META_MODE; } } } else if (attribute->global or attribute->global_early or attribute->global_late or is_global_mode()) { // Global operation: call once if ((attribute->global and do_global_mode) or (attribute->global_early and do_early_global_mode) or (attribute->global_late and do_late_global_mode) or is_global_mode()) { if (attribute->loop_local) { BEGIN_GLOBAL_MODE(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { BEGIN_LOCAL_MAP_LOOP(cctkGH, CCTK_GF) { BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Global time local mode", function, attribute, data, user_timer); } END_LOCAL_COMPONENT_LOOP; } END_LOCAL_MAP_LOOP; if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } END_REFLEVEL_LOOP; } END_GLOBAL_MODE; } else if (attribute->loop_singlemap) { BEGIN_GLOBAL_MODE(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { BEGIN_MAP_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Global time singlemap mode", function, attribute, data, user_timer); } END_MAP_LOOP; if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } END_REFLEVEL_LOOP; } END_GLOBAL_MODE; } else if (attribute->loop_level) { BEGIN_GLOBAL_MODE(cctkGH) { BEGIN_REFLEVEL_LOOP(cctkGH) { CallScheduledFunction ("Global time level mode", function, attribute, data, user_timer); if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } END_REFLEVEL_LOOP; } END_GLOBAL_MODE; } else { BEGIN_GLOBAL_MODE(cctkGH) { CallScheduledFunction ("Global mode", function, attribute, data, user_timer); if (not sync_groups.empty()) { BEGIN_REFLEVEL_LOOP(cctkGH) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } END_REFLEVEL_LOOP; } } END_GLOBAL_MODE; } } } else if (attribute->level) { // Level operation: call once per refinement level if (attribute->loop_local) { BEGIN_LOCAL_MAP_LOOP(cctkGH, CCTK_GF) { BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Level time local mode", function, attribute, data, user_timer); } END_LOCAL_COMPONENT_LOOP; } END_LOCAL_MAP_LOOP; } else if (attribute->loop_singlemap) { BEGIN_MAP_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Level time singlemap mode", function, attribute, data, user_timer); } END_MAP_LOOP; } else { CallScheduledFunction ("Level mode", function, attribute, data, user_timer); } if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } else if (attribute->singlemap) { // Single map operation: call once per refinement level and map if (attribute->loop_local) { BEGIN_LOCAL_MAP_LOOP(cctkGH, CCTK_GF) { BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Singlemap time local mode", function, attribute, data, user_timer); } END_LOCAL_COMPONENT_LOOP; } END_LOCAL_MAP_LOOP; } else { BEGIN_MAP_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Singlemap mode", function, attribute, data, user_timer); } END_MAP_LOOP; } if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } else { // Local operation: call once per component BEGIN_LOCAL_MAP_LOOP(cctkGH, CCTK_GF) { BEGIN_LOCAL_COMPONENT_LOOP(cctkGH, CCTK_GF) { CallScheduledFunction ("Local mode", function, attribute, data, user_timer); } END_LOCAL_COMPONENT_LOOP; } END_LOCAL_MAP_LOOP; if (not sync_groups.empty()) { SyncGroupsInScheduleBlock (attribute, cctkGH, sync_groups, sync_timer); } } if (schedule_barriers) { // Create an ID that is almost unique for this scheduled // function call stringstream buf; buf << cctkGH->cctk_iteration << "\n"; buf << attribute->meta << attribute->meta_early << attribute->meta_late << attribute->global << attribute->global_early << attribute->global_late << attribute->level << attribute->singlemap << attribute->local << "\n"; buf << attribute->where << "\n"; buf << attribute->thorn << "\n"; buf << attribute->routine << "\n"; string const str = buf.str(); int const id = adler32(str.c_str(), str.length()); static Timers::Timer barrier_timer ("barrier"); barrier_timer.start(); Carpet::NamedBarrier (NULL, id, "Carpet::CallFunction"); barrier_timer.stop(); } total_timer.stop(); // The return value indicates whether the grid functions have been // synchronised. // 0: let the flesh do the synchronisation // 1: we did the synchronisation return 1; } void CallScheduledFunction (char const * restrict const time_and_mode, void * const function, cFunctionData * const attribute, void * const data, Timers::Timer & user_timer) { cGH const * const cctkGH = static_cast (data); Checkpoint ("%s call at %s to %s::%s", time_and_mode, attribute->where, attribute->thorn, attribute->routine); int const skip = CallBeforeRoutines (cctkGH, function, attribute, data); if (not skip) { Timers::Timer timer(attribute->routine); // Save the time step size CCTK_REAL const saved_cctk_delta_time = cctkGH->cctk_delta_time; user_timer.start(); #ifdef REQUIREMENTS_HH Requirements::BeforeRoutine (attribute, cctkGH->cctk_iteration, reflevel, map, timelevel, timelevel_offset); #endif timer.start(); if (CCTK_IsFunctionAliased("Accelerator_PreCallFunction")) { Timers::Timer pre_timer("PreCall"); pre_timer.start(); Accelerator_PreCallFunction(cctkGH, attribute); pre_timer.stop(); } int const res = CCTK_CallFunction (function, attribute, data); assert (res==0); if (CCTK_IsFunctionAliased("Accelerator_PostCallFunction")) { Timers::Timer post_timer("PostCall"); post_timer.start(); Accelerator_PostCallFunction(cctkGH, attribute); post_timer.stop(); } timer.stop(); #ifdef REQUIREMENTS_HH Requirements::AfterRoutine (attribute, cctkGH->cctk_iteration, reflevel, map, timelevel, timelevel_offset); #endif user_timer.stop(); // Manage the time step size. If the time step size changes // during initialisation, assume it is thorn Time, and update // the time hierarchy. If it changes during evolution, assume it // is MoL, and do nothing. if (cctkGH->cctk_iteration == 0 and cctkGH->cctk_delta_time != saved_cctk_delta_time) { // The user changed cctk_delta_time during initialisation -- // update our internals and the time hierarchy bool const is_global = attribute->meta or attribute->meta_early or attribute->meta_late or attribute->global or attribute->global_early or attribute->global_late; delta_time = cctkGH->cctk_delta_time / mglevelfact * (is_global ? 1.0 : timereflevelfact); for (int ml=0; mlset_delta(ml,rl,dt); CCTK_REAL const t0 = tt->get_time(ml,rl,0); // Update the times of the past timelevels for (int tl=1; tlset_time(ml,rl,tl,t); } } } } } CallAfterRoutines (cctkGH, function, attribute, data); } void SyncGroupsInScheduleBlock (cFunctionData* attribute, cGH* cctkGH, vector const & sync_groups, Timers::Timer & sync_timer) { DECLARE_CCTK_PARAMETERS; if (sync_barriers) { // Create an ID that is almost unique for this scheduled // function call stringstream buf; buf << cctkGH->cctk_iteration << "\n"; buf << attribute->meta << attribute->meta_early << attribute->meta_late << attribute->global << attribute->global_early << attribute->global_late << attribute->level << attribute->singlemap << attribute->local << "\n"; buf << attribute->where << "\n"; buf << attribute->thorn << "\n"; buf << attribute->routine << " sync\n"; string const str = buf.str(); int const id = adler32(str.c_str(), str.length()); static Timers::Timer barrier_timer ("sync_barrier"); barrier_timer.start(); Carpet::NamedBarrier (NULL, id, "Carpet::Sync"); barrier_timer.stop(); } sync_timer.start(); SyncProlongateGroups (cctkGH, sync_groups, attribute); sync_timer.stop(); } } // namespace Carpet