aboutsummaryrefslogtreecommitdiff
path: root/src/Panda/Simple_IO.C
blob: 54998eb7b78067ec56ee1118439244cf6828bb04 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
#include "definitions.h"
#include "MPIFS.h"
#include "Chunk.h"
#include "App_Info.h"
#include "Simple_IO.h"
#include "Array.h"
#include "message.h"

#include "external/FlexIO/src/Arch.h"
#include "external/FlexIO/src/IOProtos.h"

extern MPIFS* MPIFS_global_obj;
extern int SUBCHUNK_SIZE;

extern "C" {
  int IOreserveChunk(IOFile,int,int,int*);
  int IOwriteStream(IOFile,void*,int);
  int IOreadStream(IOFile,void*,int);
  int IOwriteAttribute(IOFile, char *, int, int, void *);
  int IOreadInfo(IOFile,int *,int *,int *,int);
  int IOreadAttributeInfo(IOFile, char *,int *, int *);
  int IOreadAttribute(IOFile,int,void*);
}

/* This constructor is needed by the compute node to create a dummy object.
 * The dummy object is needed so that the compute node can execute the 
 * specialized compute node io loop
 */
Simple_IO::Simple_IO()
{
  dummy_ = YES;
  schema_string_ = current_schema_ptr_ = NULL;
  current_array_ =NULL;
  current_chunk_ = NULL;
  num_io_nodes_ = -1;
  my_io_rank_ = -1;
  compute_app_num_ = -1;
  app_info_ = NULL;
  part_time_io_ = NO;
  compute_node_array_ =NULL;
  overlap_chunk_ids_ = dest_ids_ = NULL;
  MPI_Comm_rank(MPI_COMM_WORLD, &world_rank_);
  schema_requests_ = NULL;
  requests_ =NULL;
  statuses_ =NULL;
  datatypes_ = NULL;
  schema_bufs_ = NULL;
  data_ptrs_ = NULL;
  overlap_base_ = overlap_size_ =overlap_stride_ =NULL;
  mem_buf_ = NULL;
}

Simple_IO::Simple_IO(int *schema_string, int schema_size, int world_rank,
		     int comp_app_num, int comp_app_size , App_Info *app_info,
		     IOFile fp) 
{
  int schema_buf_size;

  dummy_                = NO;
  schema_string_ 	= schema_string;
  schema_size_       	= schema_size;
  current_schema_ptr_	= schema_string;
  num_io_nodes_ 	= MPIFS_global_obj->app_size(IO_NODE);
  my_io_rank_		= MPIFS_global_obj->my_rank(IO_NODE);
  compute_app_num_	= comp_app_num;
  app_info_              = app_info;
  world_rank_            = world_rank;

  num_overlaps_         = 0;
  max_overlaps_         = comp_app_size;
  overlap_chunk_ids_    = (int *) malloc(sizeof(int)*max_overlaps_);
  dest_ids_             = (int *) malloc(sizeof(int)*max_overlaps_);
  schema_bufs_          = (int **) malloc(sizeof(int *) *max_overlaps_);
  requests_             = (MPI_Request*)malloc(sizeof(MPI_Request)*max_overlaps_);
  schema_requests_      = (MPI_Request*)malloc(sizeof(MPI_Request)*max_overlaps_);
  statuses_             = (MPI_Status*) malloc(sizeof(MPI_Status)*max_overlaps_);
  datatypes_            = (MPI_Datatype*)malloc(sizeof(MPI_Datatype)*max_overlaps_);
  max_rank_             = 10;
  overlap_base_         = (int *) malloc(sizeof(int)*max_rank_);
  overlap_stride_       = (int *) malloc(sizeof(int)*max_rank_);
  overlap_size_         = (int *) malloc(sizeof(int)*max_rank_);
  data_ptrs_            = (char **) malloc(sizeof(char*)*max_overlaps_);
  part_time_io_		= NO;
  compute_node_array_   = NULL;
  mem_buf_size_ 	= MPIFS_global_obj->mem_buf_size();
  mem_buf_              = MPIFS_global_obj->mem_buf();

  schema_buf_size = 6+ max_rank_*3;  
  for(int i=0; i < max_overlaps_; i++){
    data_ptrs_[i] = NULL;
    schema_bufs_[i] = (int *) malloc(sizeof(int)*schema_buf_size);
  }
  
  current_array_ = new Array(&schema_string);

  current_chunk_	= NULL;
  num_of_chunks_	= 0;
  num_of_subchunks_	= 0;
  current_chunk_id_ 	= -1;
  current_subchunk_id_	= -1;
  file_ptr_		= NULL;
  schema_file_ptr_	= NULL;
  file_ptr_ = fp;
}

Simple_IO::~Simple_IO()
{
  if (dummy_){
  } else {

    /* This is the object created for the I/O nodes */

    if (current_array_) delete current_array_;
    if (schema_string_) free(schema_string_);
    if (overlap_chunk_ids_) free(overlap_chunk_ids_);
    if (dest_ids_) free(dest_ids_);
    if (requests_) free(requests_);
    if (schema_requests_) free(schema_requests_);
    if (statuses_) free(statuses_);
    if (datatypes_) free(datatypes_);
    if (overlap_base_) free(overlap_base_);
    if (overlap_size_) free(overlap_size_);
    if (overlap_stride_) free(overlap_stride_);
    
    if (schema_bufs_){
      for(int i=0;i < max_overlaps_; i++){
	if (schema_bufs_[i]) free(schema_bufs_[i]);
	schema_bufs_[i] = NULL;
      }
      free(schema_bufs_);
    }

    if (data_ptrs_) free(data_ptrs_);
    
    schema_bufs_ = NULL;
    data_ptrs_ = NULL;   
    overlap_base_ = overlap_size_ = overlap_stride_ = NULL;
    overlap_chunk_ids_ = dest_ids_ = NULL;
    requests_ = NULL;
    schema_requests_ = NULL;
    statuses_ = NULL;
    datatypes_ =NULL;
    schema_string_ = NULL;	
  }
}

void Simple_IO::realloc_buffers(int new_size)
{
  int schema_buf_size = 6+max_rank_*3;


  overlap_chunk_ids_=(int *) realloc(overlap_chunk_ids_, new_size*sizeof(int));
  schema_bufs_ = (int **) realloc(schema_bufs_, new_size*sizeof(int*));
  dest_ids_ = (int *) realloc(overlap_chunk_ids_, new_size*sizeof(int));
  requests_ = (MPI_Request*)realloc(requests_, new_size*sizeof(MPI_Request));
  schema_requests_ = (MPI_Request*)realloc(schema_requests_, 
					   new_size*sizeof(MPI_Request));
  statuses_ = (MPI_Status*)realloc(statuses_, new_size*sizeof(MPI_Status));
  datatypes_ = (MPI_Datatype*)realloc(datatypes_, 
				      new_size*sizeof(MPI_Datatype));
  data_ptrs_ = (char **) realloc(data_ptrs_, new_size*sizeof(char*));
  for(int i=max_overlaps_;i<new_size;i++){
    schema_bufs_[i] = (int *)malloc(sizeof(int)*schema_buf_size);
    data_ptrs_[i] = NULL;
  }
  max_overlaps_ = new_size;
}

/* This is called only for the following cases                 *
 *  - natural chunking with user-specified subchunking         *
 *  - reorganization (with or without user-specified chunking) */
void Simple_IO::compute_chunk_overlaps(Array *array, Chunk *subchunk)
{
  int num_compute_chunks;

  if (nat_chunked_){
    num_overlaps_ = 1;
    overlap_chunk_ids_[0] = current_chunk_id_;
  }
  else{
    num_compute_chunks = array->layout(COMPUTE_NODE)->total_elements();
    if (num_compute_chunks > max_overlaps_) realloc_buffers(num_compute_chunks);
    subchunk->chunk_overlaps(array, &num_overlaps_,
			     overlap_chunk_ids_, COMPUTE_NODE);
  }

  for(int i=0; i < num_overlaps_;i++) {
    dest_ids_[i]=app_info_->world_rank(array->which_node(overlap_chunk_ids_[i],
						         COMPUTE_NODE));
}

#ifdef DEBUG
  printf("For subchunk_id %d of chunk %d\n", current_subchunk_id_, 
					     current_chunk_id_);
  printf("The overlapping compute chunk ids are \n");
  for(int k =0; k < num_overlaps_; k++) printf("%d ", overlap_chunk_ids_[k]);
  printf("\n");
#endif
}


/* This is called only for the following cases                 *
 *  - natural chunking with user-specified subchunking         *
 *  - reorganization (with or without user-specified chunking) */
void Simple_IO::compute_schemas(Array *array, Chunk *subchunk , 
				Chunk *compute_chunk)
{
  if (nat_chunked_ && !contiguous_ && !overlaped_){
    subchunk->copy_base_size_stride(overlap_base_, overlap_size_, 
				    overlap_stride_);
    send_schema_message(0);
    make_datatype(subchunk, 0);
  }
  else if (!nat_chunked_) {
    for (int i=0; i< num_overlaps_; i++){
      compute_chunk->init(array, overlap_chunk_ids_[i], COMPUTE_NODE, NO_ALLOC);
      subchunk->compute_overlap(compute_chunk, overlap_base_, overlap_size_, 
 				overlap_stride_);
      send_schema_message(i);
      make_datatype(subchunk, i);
    }
  } else {
    printf("Error - In Simple_IO::compute_schemas\n");
    exit(1);
  }
}


/* The chunk_id is in overlap_chunk_ids_[index], the dest is in  *
 * in dest_ids_[index]. The rank,base,stride and size info is in *
 * overlap_base, overlap_size, overlap_stride, array_rank_       */
void Simple_IO::send_schema_message(int index)
{
   int *ptr = schema_bufs_[index];
   int schema_size = 5+array_rank_*3;

   *ptr++ = overlap_chunk_ids_[index];
   *ptr++ = (int) nat_chunked_;
   *ptr++ = (int) contiguous_;
   *ptr++ = array_rank_;
   *ptr++ = op_type_;

   for(int i=0; i < array_rank_; i++) *ptr++ = overlap_base_[i];
   for(i=0; i < array_rank_; i++) *ptr++ = overlap_size_[i];
   for(i=0; i < array_rank_; i++) *ptr++ = overlap_stride_[i];

   if (part_time_io_ && (dest_ids_[index] == world_rank_))
     /* No need to send the  message */
     schema_requests_[index] = MPI_REQUEST_NULL;
   else 
     nb_send_message((void *)schema_bufs_[index], schema_size, MPI_INT,
		     dest_ids_[index], index*10+CHUNK_SCHEMA, MPI_COMM_WORLD,
		     &schema_requests_[index]);
}

/* The overlap base, size, stride are in overlap_base, overlap_size, *
 * and overlap_stride                                                */
void Simple_IO::make_datatype(Chunk *subchunk, int index)
{
  void *ptr;
  subchunk->make_datatype(overlap_base_, overlap_size_, overlap_stride_, 
			  &ptr, &datatypes_[index]);
  data_ptrs_[index] = (char *) ptr;
}

/* Again this function is called only for the following cases  *
 * - natural chunking with user-specified subchunking          *
 * - re-organization with/without user-specified chunking      *
 * The case of natural chunking (with no user-specified        *
 * subchunking) is handled seperately                          */

void Simple_IO::receive_data(Chunk *subchunk, int index, int &array_bytes_to_go)
{

  if (part_time_io_ && (dest_ids_[index] == world_rank_)){
    /* Perform a mem copy of the required chunk */
    copy_data(subchunk, index, NO, array_bytes_to_go);
    requests_[index] = MPI_REQUEST_NULL;
  } else 
    nb_receive_message((void *)data_ptrs_[index], 1, datatypes_[index],
		       dest_ids_[index], index*10+CHUNK_DATA_TO_IO, 
		       MPI_COMM_WORLD, &requests_[index]);
}

/* Again this function is called only for the following cases  *
 * - natural chunking with user-specified subchunking          *
 * - re-organization with/without user-specified chunking      *
 * The case of natural chunking (with no user-specified        *
 * subchunking) is handled seperately                          */
void Simple_IO::send_data(Chunk *subchunk, int index, int &array_bytes_to_go)
{
  if (part_time_io_ && (dest_ids_[index] == world_rank_)){
    /* Perform a memory copy of the required chunk */
    copy_data(subchunk, index, YES, array_bytes_to_go);
    requests_[index] =MPI_REQUEST_NULL;
  } else {
    /* Send the required datatype using a non-blocking send */
    nb_send_message((void *)data_ptrs_[index], 1, datatypes_[index],
		    dest_ids_[index], index*10+CHUNK_DATA_FROM_IO, 
		    MPI_COMM_WORLD, &requests_[index]);
  }
}
  
void Simple_IO::read_data(Chunk *subchunk)
{
  int size;
  size = subchunk->total_size_in_bytes();
  read_data((char *)(subchunk->data_ptr()), size, subchunk->element_size());
}     

void Simple_IO::read_data(char *buf, int size, int esize)
{
  int n,bytes_to_go=size,buf_size;
  char *tmp_buf = buf;

  while(bytes_to_go > 0){
    buf_size = min(bytes_to_go, SUBCHUNK_SIZE);
    n = IOreadStream(file_ptr_, (void *)tmp_buf, buf_size/esize);
    if (n != buf_size){
      printf("Error reading data - write only %d instead of %d bytes\n", 
	      n, buf_size);
//    exit(1);
    } 
    bytes_to_go -= buf_size;
    tmp_buf += buf_size;
  }
}


void Simple_IO::write_data(char *buf, int size, int esize)
{
  int n, bytes_to_go = size, buf_size;
  char *tmp_buf = buf;

  while(bytes_to_go > 0){
    buf_size = min(bytes_to_go, SUBCHUNK_SIZE);
    n = IOwriteStream(file_ptr_, (void *)tmp_buf, buf_size/esize);
    if (n != buf_size){
      printf("Error writing data - write only %d instead of %d bytes\n",
	     n, buf_size);
      exit(1);
    } 
    tmp_buf += buf_size;
    bytes_to_go -= buf_size;
  }
}

void Simple_IO::write_data(Chunk* subchunk)
{
  int size;
  size = subchunk->total_size_in_bytes();
  write_data((char *)(subchunk->data_ptr()), size, subchunk->element_size());
}     

void Simple_IO::free_datatypes()
{
 for(int i=0; i <num_overlaps_; i++) MPI_Type_free(&datatypes_[i]);
}

void Simple_IO::send_data_to_compute_nodes(Chunk *subchunk, 
					   int &array_bytes_to_go)
{
  for(int i=0; i< num_overlaps_; i++)
    send_data(subchunk, i, array_bytes_to_go);
}

void Simple_IO::receive_data_from_compute_nodes(Chunk *subchunk, 
						int &array_bytes_to_go)
{
  for (int i=0; i< num_overlaps_; i++) 
    receive_data(subchunk, i, array_bytes_to_go);
}

void Simple_IO::wait_for_completion(int &array_bytes_to_go,
				    Array *compute_array)
{
  int flag=0;

  if (part_time_io_){
    /* This is to avoid deadlocks */
    while (!flag){
      MPI_Testall(num_overlaps_, requests_, &flag, statuses_);
      if (array_bytes_to_go > 0) 
	process_compute_message(array_bytes_to_go, compute_array);
    }
  } else {
    MPI_Waitall(num_overlaps_, requests_, statuses_);
  }
  /* Free the schema request objects - Do we need this*/
  MPI_Waitall(num_overlaps_, schema_requests_, statuses_);
}

/* For part-io nodes, get the data using memory copy if the *
 * data resides on the same node.                           */
void Simple_IO::copy_data(Chunk *subchunk, int index, Boolean flag,
			  int &array_bytes_to_go)
{
  void *comp_data_ptr;
  MPI_Datatype comp_datatype;
  int position=0, buf_size;
  void *buf=NULL;
  int *schema = schema_bufs_[index];
  int comp_chunk_id =  schema[0];
  int comp_array_rank = schema[3];
  int *base = &schema[5];
  int *size = &schema[5+comp_array_rank*1];
  int *stride = &schema[5+comp_array_rank*2];
  int bytes_copied = num_elements(comp_array_rank, size)*
                     subchunk->element_size();
  Array *comp_array = compute_node_array_;
  Chunk *comp_chunk = comp_array->find_chunk(comp_chunk_id);
  comp_chunk->make_datatype(base, size,stride, &comp_data_ptr, 
			    &comp_datatype);
  if (array_bytes_to_go > 0) array_bytes_to_go -= bytes_copied;

  if (flag){
    MPI_Pack_size(1, datatypes_[index], MPI_COMM_WORLD, &buf_size);
    buf = (void *) malloc(buf_size);
    MPI_Pack(data_ptrs_[index], 1, datatypes_[index], buf, buf_size, 
	     &position, MPI_COMM_WORLD);
    position =0;
    MPI_Unpack(buf, buf_size, &position, comp_data_ptr, 1, comp_datatype, 
	       MPI_COMM_WORLD);
    free(buf);
  } else {
    MPI_Pack_size(1, comp_datatype, MPI_COMM_WORLD, &buf_size);
    buf = (void *) malloc(buf_size);
    MPI_Pack(comp_data_ptr, 1, comp_datatype, buf, buf_size, 
	     &position, MPI_COMM_WORLD);
    position = 0;
    MPI_Unpack(buf, buf_size, &position, data_ptrs_[index], 1,
	       datatypes_[index], MPI_COMM_WORLD);
    free(buf);
  }
  MPI_Type_free(&comp_datatype);
}

/* For nat chunking with no user defined subchunking, read/write
 * data directly from compute chunk (i.e if it is on same node) */    
void Simple_IO::direct_io(int chunk_id, Boolean flag, int &array_bytes_to_go)
{
  Array *comp_array = compute_node_array_;
  Chunk *comp_chunk = comp_array->find_chunk(chunk_id);
  if (flag) read_data(comp_chunk);
  else write_data(comp_chunk);
  if (array_bytes_to_go > 0) 
    array_bytes_to_go -= comp_chunk->total_size_in_bytes();
}

void Simple_IO::realloc_schema_bufs(int new_size)
{
  int schema_buf_size = sizeof(int)*(6+new_size*3);
  
  max_rank_ = new_size;
  overlap_base_ = (int *) realloc(overlap_base_, max_rank_*sizeof(int));
  overlap_size_ = (int *) realloc(overlap_size_, max_rank_*sizeof(int));
  overlap_stride_ = (int *) realloc(overlap_stride_, max_rank_*sizeof(int));
  for(int i=0; i < max_overlaps_; i++){
    schema_bufs_[i] = (int *) realloc(schema_bufs_[i], schema_buf_size);
  }
}
  
void Simple_IO::realloc_mem_bufs(int new_size)
{
  mem_buf_size_ = new_size;
  mem_buf_ = (char *) realloc(mem_buf_, sizeof(char)*mem_buf_size_);
  MPIFS_global_obj->set_mem_buf_size(new_size);
  MPIFS_global_obj->set_mem_buf(mem_buf_);
}

void Simple_IO::start_to_finish(Boolean part_time, Array *compute_array)
{
  int make_subchunks, bytes_to_go;
  int array_bytes_to_go,*ptr;
  Boolean read_op;
  Chunk *chunk=NULL, *subchunk=NULL, *compute_chunk=NULL, *tmp_chunk;

  op_type_ 	 = current_array_->op_type();
  if ((op_type_ == RESTART)||(op_type_ == GENERAL_READ)||
      (op_type_ == READ_TIMESTEP))
    read_op = YES;
  else
    read_op = NO;

  part_time_io_ = part_time;
  compute_node_array_ = compute_array;

  if (read_op) {
    int numbertype, rank, index, datatype, length;
    int *dims = (int *)malloc(sizeof(int) * 10);
    IOreadInfo(file_ptr_, &numbertype, &rank, dims, 10);
    int *size = (int *)malloc(sizeof(int) * rank);

    index = IOreadAttributeInfo(file_ptr_, "global_size", &datatype, &length);
    if (index >=0 ) { // the attribute exists
      IOreadAttribute(file_ptr_, index, size);
      current_array_->init(rank, numbertype, size, IO_NODE); 
    } else { printf("Error: no attribute, global_size\n"); exit(0); }
    free(dims);

printf("%d: read rank %d, numbertype %d, size (%d %d %d)\n", world_rank_,
	rank, numbertype, size[0], size[1], size[2]);

    int schema_size = 2 + rank;
    int *schema = (int *)malloc(sizeof(int) * schema_size);
    if (MPIFS_global_obj->am_master_io_node()) {
      schema[0] = rank; schema[1] = numbertype; 
      for (int i=0; i<rank; i++) schema[2+i] = size[i];
      send_message((void *)schema, schema_size, MPI_INT, 
		   app_info_->get_master(),
		   ARRAYGROUP_SCHEMA, MPI_COMM_WORLD);
    }
    if (part_time_io_) {
      MPI_Status status;
      receive_message(schema, schema_size, MPI_INT, MPI_ANY_SOURCE, 
		      ARRAYGROUP_SCHEMA, MPI_COMM_WORLD, &status);
      MPIFS_global_obj->Broadcast(COMPUTE_NODE, (void *)schema,
                                  schema_size, MPI_INT, ARRAYGROUP_SCHEMA);

      compute_array->init(rank, numbertype, size, COMPUTE_NODE);
    }
    free(schema);
  }
 
  if (part_time_io_) array_bytes_to_go = compute_node_array_->array_info();
 
  /* To reduce costs associated with object creation and deletion, we *
   * will create a dummy chunk,subchunk and compute chunk object and  *
   * re-initialize them whenever necessary.                           */
  tmp_chunk = chunk = new Chunk();
  current_chunk_ = chunk;
  subchunk = new Chunk();
  compute_chunk = new Chunk();

  make_subchunks = -1;

  nat_chunked_   = current_array_->nat_chunked();
  sub_chunked_   = current_array_->sub_chunked();
  overlaped_     = current_array_->overlaped();
  if (overlaped_) { contiguous_ = NO; nat_chunked_ = NO; }
  else {  
    if (nat_chunked_ && !sub_chunked_)
      contiguous_ = YES; /* No need to use derived datatypes */
    else contiguous_ = NO;  /* Have to use derived datatypes */
  }
    
  array_rank_ = current_array_->rank();
  if (array_rank_ > max_rank_) realloc_schema_bufs(array_rank_);

  if (read_op) current_array_->read_schema_file(file_ptr_);

  num_of_chunks_ = current_array_->layout(IO_NODE)->total_elements();
  current_chunk_id_ = current_array_->get_next_index(chunk, -1, my_io_rank_,
						     num_io_nodes_,
						     num_of_chunks_);

#ifdef DEBUG
    printf("%d: current_chunk_id_=%d my_io_rank=%d num_io_nodes=%d\n",
	   world_rank_, current_chunk_id_, my_io_rank_, num_io_nodes_);
#endif
  if (contiguous_){
    /* Natural chunked and no user-specified subchunking. Therefore we don't 
     *  need to used mpi-derived datatypes.  */
      
    while (current_chunk_id_ < num_of_chunks_) {
      if (!read_op) { 
        int *tmp_size = (int *)malloc(sizeof(int) * array_rank_);
        for (int cnt = 0; cnt < array_rank_; cnt++)
          tmp_size[cnt] = chunk->size()[array_rank_ - cnt - 1];
        IOreserveChunk(file_ptr_, current_array_->ieee_size(), 
		       array_rank_, tmp_size);
        //printf("##### called IOreserveChunk for n.c. %d %d %d %d %d\n", current_array_->ieee_size(), array_rank_, tmp_size[0], tmp_size[1], tmp_size[2]);
 
        free(tmp_size);
        if (num_of_chunks_ > 1) {
	  IOwriteAttribute(file_ptr_,"chunk_origin", INT32, 3, chunk->base());
	  IOwriteAttribute(file_ptr_, "chunk_size", INT32, 3, chunk->size());
	}
      } 

      /* for part-time io case, if chunk resides on same node, perform the *
       * read/write operation directly.                                    */
      num_overlaps_ = 1;
      overlap_chunk_ids_[0] = current_chunk_id_;
      dest_ids_[0] = app_info_->world_rank(current_array_->which_node(
					   current_chunk_id_, COMPUTE_NODE));

      if (part_time_io_ && (world_rank_ == dest_ids_[0])){
	direct_io(current_chunk_id_, read_op, array_bytes_to_go);
      } else {
	bytes_to_go = chunk->total_size_in_bytes();
	chunk->set_data_ptr(mem_buf_);

	/* Make the schema request */
	ptr = schema_bufs_[0];
	*ptr++ = current_chunk_id_;
	*ptr++ = (int)nat_chunked_;
	*ptr++ = (int)contiguous_;
	*ptr++ = op_type_;
	*ptr++ = 0;  /* This is the offset */
	*ptr++ = 0; /* Size of the data */

	ptr = schema_bufs_[0];
	while(bytes_to_go > 0){
	  ptr[5] = min(SUBCHUNK_SIZE, bytes_to_go);

	  nb_send_message((void *)ptr, 6, MPI_INT, dest_ids_[0],
			  CHUNK_SCHEMA, MPI_COMM_WORLD, &schema_requests_[0]);
	  if (read_op){
	    read_data(mem_buf_, ptr[5], chunk->element_size());
	    nb_send_message((void *)mem_buf_, ptr[5], MPI_CHAR, dest_ids_[0],
			    CHUNK_DATA_FROM_IO, MPI_COMM_WORLD, &requests_[0]);
	  } else
	    nb_receive_message((void *)mem_buf_, ptr[5], MPI_CHAR, 
				dest_ids_[0], CHUNK_DATA_TO_IO, 
				MPI_COMM_WORLD, &requests_[0]);
	  /* Have to watch for deadlock over here */
	  wait_for_completion(array_bytes_to_go, compute_node_array_);
	  if (!read_op) write_data(mem_buf_, ptr[5], chunk->element_size());
	  ptr[4] += ptr[5];
	  bytes_to_go -= ptr[5];
	}
	chunk->set_data_ptr(NULL);
      }
      current_chunk_id_ = current_array_->get_next_index(chunk, 
							 current_chunk_id_,
							 my_io_rank_,
							 num_io_nodes_,
							 num_of_chunks_);
    } /* End while */
  } /* End if (contiguous_) */
  else {
    /* We have no choice but to use MPI-derived datatypes */
    while(current_chunk_id_ < num_of_chunks_){
      if (!read_op) { 
        int *tmp_size = (int *)malloc(sizeof(int) * array_rank_);
        for (int cnt = 0; cnt < array_rank_; cnt++)
          tmp_size[cnt] = chunk->size()[array_rank_ - cnt - 1];
        IOreserveChunk(file_ptr_, current_array_->ieee_size(), 
		array_rank_, tmp_size);
        //printf("##### called IOreserveChunk for r.o. %d %d %d %d %d\n", current_array_->ieee_size(), array_rank_, tmp_size[0], tmp_size[1], tmp_size[2]);

        free(tmp_size);
        if (num_of_chunks_ > 1) {
          IOwriteAttribute(file_ptr_,"chunk_origin", INT32, 3, chunk->base());
	  IOwriteAttribute(file_ptr_, "chunk_size", INT32, 3, chunk->size());
        }
      } 

      /* If the array is not subchunked, then subchunk the array into *
       * SUBCHUNK_SIZE chunks. This is to reduce the size of the      *
       * messages and the memory requirements. The current version makes a   *
       * dumb assumption, that if the user specifies the subchunks,   *
       * then the size of those subchunks is less than SUBCHUNK_SIZE. *
       * It's a dumb assumption and needs to be fixed.                */

      if (!sub_chunked_ && (make_subchunks == -1)){
        current_array_->make_sub_chunks(chunk);
	make_subchunks = 1;
      }
      num_of_subchunks_ =current_array_->layout(SUB_CHUNK)->total_elements();

      for (current_subchunk_id_=0; current_subchunk_id_ < num_of_subchunks_;
	   current_subchunk_id_++){
	subchunk->init(chunk, current_subchunk_id_, NO_ALLOC);
	bytes_to_go = subchunk->total_size_in_bytes();

	if (bytes_to_go > mem_buf_size_) realloc_mem_bufs(bytes_to_go);
	subchunk->set_data_ptr(mem_buf_);

	compute_chunk_overlaps(current_array_, subchunk);
	compute_schemas(current_array_, subchunk, compute_chunk);

        if (read_op){
	  read_data(subchunk);
	  send_data_to_compute_nodes(subchunk, array_bytes_to_go);
	} else receive_data_from_compute_nodes(subchunk, array_bytes_to_go);
	wait_for_completion(array_bytes_to_go, compute_node_array_);
	if (!read_op) write_data(subchunk);

	free_datatypes();
	subchunk->set_data_ptr(NULL);
      }
      current_chunk_id_ = current_array_->get_next_index(chunk, 
							 current_chunk_id_,
							 my_io_rank_,
							 num_io_nodes_,
							 num_of_chunks_);
    } /* End while loop */
  } /* End if else */

#ifdef DEBUG
  printf("%d:Finished the I/O\n", world_rank_);
#endif
  if (part_time_io_){
    /* Since the I/O side is finished jump into the compute loop */
    while (array_bytes_to_go > 0)
      process_compute_message(array_bytes_to_go, compute_node_array_);
#ifdef DEBUG
    printf("%d:Finished the compute side of the part-time io\n", world_rank_);
#endif
  }

  /* Delete chunk, subchunk, compute_chunk */
  if (tmp_chunk) delete tmp_chunk;
  if (subchunk) delete subchunk;
  if (compute_chunk) delete compute_chunk;
  chunk=subchunk=compute_chunk=NULL;
}

/* This function should not access any of the instance variables of
 * the Simple_IO object without setting them first
 */
void Simple_IO::compute_node_io_loop(Array *array)
{
  int op_type = array->op_type();
  if ((op_type == RESTART) || (op_type == GENERAL_READ) ||
      (op_type == READ_TIMESTEP)) {
    MPI_Status status;
    int *schema, schema_size;

    MPI_Probe(MPI_ANY_SOURCE, ARRAYGROUP_SCHEMA, MPI_COMM_WORLD, &status);
    mpi_get_count(&status, MPI_INT, &schema_size);
    schema = (int *)malloc(sizeof(int) * schema_size);
    receive_message((void *)schema, schema_size, MPI_INT, status.MPI_SOURCE,
		    ARRAYGROUP_SCHEMA, MPI_COMM_WORLD, &status);
    MPIFS_global_obj->Broadcast(COMPUTE_NODE, (void *)schema,
                                schema_size, MPI_INT, ARRAYGROUP_SCHEMA);

    int *size = (int *)malloc(sizeof(int) * schema[0]);
    for (int i=0; i<schema[0]; i++) size[i] = schema[2+i];
printf("%d: read rank %d, numbertype %d, size (%d %d %d)\n", world_rank_,
	schema[0], schema[1], size[0], size[1], size[2]);
    array->init(schema[0], schema[1], size, COMPUTE_NODE);
    free(schema);
  }

  int array_bytes_to_go = array->array_info();
  while (array_bytes_to_go > 0) 
    process_compute_message(array_bytes_to_go, array);
}

void Simple_IO::process_compute_message(int &arrays_bytes_to_go, 
				        Array *array)
{
  int msg_code, msg_tag, msg_src;
  MPI_Status status;
  int data_size;

  any_new_message(&msg_code, &msg_src, &msg_tag, &status);
  switch(msg_code){
  case CHUNK_SCHEMA:
    /* Do something about it */
    process_chunk_schema_request(msg_src,msg_tag, arrays_bytes_to_go, 
				 &status, array);
    break;
    
  case CHUNK_DATA_FROM_IO:
    MPI_Get_count(&status, MPI_CHAR, &data_size);
    printf("Received chunk_data before chunk schema from %d of size %d\n",
	   msg_src, data_size);
    MPI_Probe(msg_src, (msg_tag/10)*10+CHUNK_SCHEMA, MPI_COMM_WORLD, &status);
    printf("Received the corressponding chunk schema message\n");
    process_chunk_schema_request(msg_src, (msg_tag/10)*10+CHUNK_SCHEMA,
				 arrays_bytes_to_go,
				 &status, array);
    break;

  case NO_MESSAGE:
    /* Do nothing */
    break;
  default:
    /* This message is not for me */
    printf("In process compute message - unknown code %d\n", msg_code);
    break;
  }
}

void Simple_IO::process_chunk_schema_request(int msg_src, int msg_tag, 
					     int &array_bytes_to_go,
					     MPI_Status *status, Array *array)
{
  int *schema_buf, schema_size;
  int chunk_id, op_type, array_rank, *base, *size, *stride, *ptr;
  int data_size, elt_size, offset;
  Boolean contiguous;
  MPI_Datatype datatype;
  Chunk *chunk;
  void *data_ptr;

  MPI_Get_count(status, MPI_INT, &schema_size);
  schema_buf = (int *) malloc(sizeof(int)*schema_size);
  receive_message((void *)schema_buf, schema_size, MPI_INT, msg_src, 
		  msg_tag, MPI_COMM_WORLD, status);

  ptr = schema_buf;
  chunk_id = *ptr++;
  ptr++;
  contiguous = (Boolean) *ptr++;
  chunk = array->find_chunk(chunk_id);

  if (contiguous){
    op_type = *ptr++;
    offset =  *ptr++;
    data_size = *ptr++;
    data_ptr = chunk->data_ptr();
    data_ptr = (char *)((char *) data_ptr + offset);
      
    if ((op_type == RESTART) || (op_type == READ_TIMESTEP) ||
	(op_type == GENERAL_READ))
      receive_message((void *) data_ptr, 
		      data_size,
		      MPI_CHAR, msg_src, 
		      (msg_tag/10*10)+CHUNK_DATA_FROM_IO,
		      MPI_COMM_WORLD, status);
      else
	send_message((void *)data_ptr, data_size,MPI_CHAR, msg_src,
		     (msg_tag/10)*10+CHUNK_DATA_TO_IO,
		     MPI_COMM_WORLD);
    
  }
  else{
    array_rank = *ptr++;
    op_type = *ptr++;
    base = &ptr[0];
    size = &ptr[array_rank*1];
    stride = &ptr[array_rank*2];
    elt_size = chunk->element_size();
    data_size =  num_elements(array_rank, size)*elt_size;
      
    chunk->make_datatype(base,size,stride, &data_ptr, &datatype);
    if ((op_type == RESTART) || (op_type == READ_TIMESTEP) ||
	(op_type == GENERAL_READ))
      receive_message(data_ptr, 1, datatype,msg_src, 
		      (msg_tag/10)*10+CHUNK_DATA_FROM_IO,
		      MPI_COMM_WORLD, status);
    else
      send_message(data_ptr, 1, datatype, msg_src,
		   (msg_tag/10)*10+CHUNK_DATA_TO_IO,
		   MPI_COMM_WORLD); 
      MPI_Type_free(&datatype);
  }
  
  array_bytes_to_go -= data_size;
  free(schema_buf);
}