5 #include "../shared/lapack_symbs.h" 6 #include "../tensor/algstrct.h" 7 #include "../shared/util.h" 8 #include "../shared/memcontrol.h" 9 #include "../shared/offload.h" 20 if (this->
pr == other.
pr)
21 return this->
pc < other.
pc;
23 return this->
pr < other.
pr;
36 comm = MPI_COMM_WORLD;
71 this->init(
comm, order, lens, argc, argv);
78 printf(
"CTF WARNING: Creating copy of World, which is not free or useful, pass original World by reference instead if possible.\n");
98 for (
int i=0; i<(int)topovec.size(); i++){
101 delete phys_topology;
102 if (this->cdt.cm == MPI_COMM_WORLD){
104 universe_exists =
false;
129 int World::init(MPI_Comm
const global_context,
132 const char *
const * argv){
135 phys_topology = NULL;
139 return initialize(argc, argv);
142 int World::init(MPI_Comm
const global_context,
146 const char *
const * argv){
149 phys_topology =
new topology(order, dim_len, cdt, 1);
151 return initialize(argc, argv);
155 int World::initialize(
int argc,
156 const char *
const * argv){
157 char * mst_size, * stack_size, *
mem_size, * ppn;
158 if (
comm == MPI_COMM_WORLD && universe_exists){
159 delete phys_topology;
167 if (phys_topology == NULL){
186 MPI_Comm_size(MPI_COMM_WORLD, &all_np);
189 printf(
"CTF ERROR: the first CTF instance created has to be on MPI_COMM_WORLD\n");
201 char * ntd = getenv(
"OMP_NUM_THREADS");
203 omp_set_num_threads(1);
205 VPRINTF(1,
"Running with 1 thread using omp_set_num_threads(1), because OMP_NUM_THREADS is not defined\n");
208 if (
rank == 0 && ntd != NULL){
209 VPRINTF(1,
"Running with %d threads\n",omp_get_max_threads());
214 char * file_path = getenv(
"CTF_MODEL_FILE");
215 if (file_path != NULL && strcmp(file_path,
"")!=0){
216 VPRINTF(1,
"Reading model coefficients from file %s (CTF_MODEL_FILE)\n", file_path);
222 mst_size = getenv(
"CTF_MST_SIZE");
223 stack_size = getenv(
"CTF_STACK_SIZE");
224 if (mst_size == NULL && stack_size == NULL){
227 VPRINTF(1,
"Creating stack of size %ld\n",1000*(int64_t)1E6);
236 int64_t imst_size = 0 ;
237 if (mst_size != NULL)
238 imst_size = strtoull(mst_size,NULL,0);
239 if (stack_size != NULL)
240 imst_size =
MAX(imst_size,strtoull(stack_size,NULL,0));
242 printf(
"Creating stack of size %ld due to CTF_STACK_SIZE enviroment variable\n",
247 mem_size = getenv(
"CTF_MEMORY_SIZE");
248 if (mem_size != NULL){
249 int64_t imem_size = strtoull(mem_size,NULL,0);
251 VPRINTF(1,
"Memory size set to %ld by CTF_MEMORY_SIZE environment variable\n",
255 ppn = getenv(
"CTF_PPN");
258 printf(
"Assuming %d processes per node due to CTF_PPN environment variable\n",
271 if (
comm == MPI_COMM_WORLD){
272 if (!universe_exists){
273 universe_exists =
true;
310 if (!universe_exists){
313 delete pscp_universe;
void load_all_models(std::string file_name)
void set_main_args(int argc, const char *const *argv)
void mst_create(int64_t size)
initializes stack buffer
an instance of the CTF library (world) on a MPI communicator
void mem_exit(int rank)
exit instance of memory manager
void init_rng(int rank)
initialized random number generator
std::vector< topology * > peel_perm_torus(topology *phys_topology, CommData cdt)
folds specified topology and all of its permutations into all configurations of lesser dimensionality...
int rank
rank of local processor
void set_context(MPI_Comm ctxt)
std::set< grid_wrapper > scalapack_grids
index for ScaLAPACK processor grids
CTF_int::topology * phys_topology
main torus topology corresponding to the world
topology * get_phys_topo(CommData glb_comm, TOPOLOGY mach)
get dimension and torus lengths of specified topology
void offload_exit()
exit offloading, e.g. destroy cublas
void cblacs_gridexit(int contxt)
void set_memcap(double cap)
sets what fraction of the memory capacity CTF can use
~World()
frees CTF library
void mem_create()
create instance of memory manager
int64_t proc_bytes_available()
gives total memory available on this MPI process
World(int argc, char *const *argv)
creates CTF library on comm that can output profile data into a file with a name based on the main ar...
std::vector< topology * > get_generic_topovec(CommData cdt)
computes all topology configurations given undelying physical topology information ...
void offload_init()
initialize offloading, e.g. create cublas
bool operator<(grid_wrapper const &other) const
void set_mem_size(int64_t size)
sets what fraction of the memory capacity CTF can use
MPI_Comm comm
set of processors making up this world