Ejemplo n.º 1
0
void test_create_simple() {
  test_work_area_type * work_area = test_work_area_alloc("nnc-INIT");
  {
    int nx = 10;
    int ny = 10;
    int nz = 10;
    ecl_grid_type * grid0 = ecl_grid_alloc_rectangular(nx,ny,nz,1,1,1,NULL);

    ecl_grid_add_self_nnc(grid0, 0 ,nx*ny + 0, 0 );
    ecl_grid_add_self_nnc(grid0, 1 ,nx*ny + 1, 1 );
    ecl_grid_add_self_nnc(grid0, 2 ,nx*ny + 2, 2 );
    {
      ecl_nnc_geometry_type * nnc_geo = ecl_nnc_geometry_alloc( grid0 );
      test_assert_int_equal( ecl_nnc_geometry_size( nnc_geo ) , 3 );

      /*
        Create a dummy INIT file which *ony* contains a TRANNC keyword with the correct size.
      */
      {
        ecl_kw_type * trann_nnc = ecl_kw_alloc(TRANNNC_KW , ecl_nnc_geometry_size( nnc_geo ), ECL_FLOAT);
        fortio_type * f = fortio_open_writer( "TEST.INIT" , false, ECL_ENDIAN_FLIP );

        for (int i=0; i < ecl_kw_get_size( trann_nnc); i++)
          ecl_kw_iset_float( trann_nnc , i , i*1.0 );

        ecl_kw_fwrite( trann_nnc , f );
        fortio_fclose( f );
        ecl_kw_free( trann_nnc );
      }
    }
    ecl_grid_free( grid0 );
  }
  test_work_area_free( work_area );
}
Ejemplo n.º 2
0
void test_kw_io_charlength() {
  test_work_area_type * work_area = test_work_area_alloc("ecl_kw_io_charlength");
  { 
    const char * KW0 = "QWERTYUI";
    const char * KW1 = "ABCDEFGHIJTTTTTTTTTTTTTTTTTTTTTTABCDEFGHIJKLMNOP";
    ecl_kw_type * ecl_kw_out0 = ecl_kw_alloc(KW0 , 5, ECL_FLOAT);
    ecl_kw_type * ecl_kw_out1 = ecl_kw_alloc(KW1 , 5, ECL_FLOAT);
    for (int i=0; i < ecl_kw_get_size( ecl_kw_out1); i++) {
       ecl_kw_iset_float( ecl_kw_out0 , i , i*1.5 );
       ecl_kw_iset_float( ecl_kw_out1 , i , i*1.5 );
    }

    {
       fortio_type * f = fortio_open_writer( "TEST1" , false, ECL_ENDIAN_FLIP );
       test_assert_true( ecl_kw_fwrite( ecl_kw_out0, f ));
       test_assert_false(ecl_kw_fwrite( ecl_kw_out1, f ));
       fortio_fclose( f );
    }

    {
       test_assert_false( util_file_exists( "TEST1"));
    }
   
    {
       FILE * file = util_fopen("TEST2", "w");
       ecl_kw_fprintf_grdecl(ecl_kw_out1 , file);
       fclose(file);
    }
    
    {
       FILE * file = util_fopen("TEST2", "r");
       ecl_kw_type * ecl_kw_in = ecl_kw_fscanf_alloc_grdecl( file , KW1 , -1 , ECL_FLOAT);
       test_assert_string_equal(KW1, ecl_kw_get_header(ecl_kw_in) );
       test_assert_int_equal(5, ecl_kw_get_size( ecl_kw_in) );

       test_assert_double_equal(ecl_kw_iget_as_double(ecl_kw_in, 0), 0.0);
       test_assert_double_equal(ecl_kw_iget_as_double(ecl_kw_in, 4), 6.0);

       fclose(file);
    }
    
    ecl_kw_free( ecl_kw_out0 );
    ecl_kw_free( ecl_kw_out1 );
  }
  test_work_area_free( work_area );
}
Ejemplo n.º 3
0
static void ecl_init_file_fwrite_poro( fortio_type * fortio , const ecl_grid_type * ecl_grid , const ecl_kw_type * poro ) {
  {
    ecl_kw_type * porv = ecl_kw_alloc( PORV_KW , ecl_grid_get_global_size( ecl_grid ) , ECL_FLOAT_TYPE);
    int global_index;
    bool global_poro = (ecl_kw_get_size( poro ) == ecl_grid_get_global_size( ecl_grid )) ? true : false;
    for ( global_index = 0; global_index < ecl_grid_get_global_size( ecl_grid ); global_index++) {
      int active_index = ecl_grid_get_active_index1( ecl_grid , global_index );
      if (active_index >= 0) {
        int poro_index = global_poro ? global_index : active_index;
        ecl_kw_iset_float( porv , global_index , ecl_kw_iget_float( poro , poro_index ) * ecl_grid_get_cell_volume1( ecl_grid , global_index ));
      } else
        ecl_kw_iset_float( porv , global_index , 0 );
    }
    ecl_kw_fwrite( porv , fortio );
    ecl_kw_free( porv );
  }

  ecl_kw_fwrite( poro , fortio );
}
Ejemplo n.º 4
0
Archivo: ecl_file.c Proyecto: akva2/ert
void test_writable(const char * src_file ) {
  test_work_area_type * work_area = test_work_area_alloc("ecl_file_writable" );
  char * fname = util_split_alloc_filename( src_file );

  test_work_area_copy_file( work_area , src_file );
  {
    ecl_file_type * ecl_file = ecl_file_open( fname , ECL_FILE_WRITABLE);
    ecl_kw_type * swat = ecl_file_iget_named_kw( ecl_file , "SWAT" , 0 );
    ecl_kw_type * swat0 = ecl_kw_alloc_copy( swat );
    test_assert_true( ecl_kw_equal( swat , swat0 ));
    ecl_kw_iset_float( swat , 0 , 1000.0 );
    ecl_file_save_kw( ecl_file , swat );
    test_assert_true( ecl_file_writable( ecl_file ));
    ecl_file_close( ecl_file );

    ecl_file = ecl_file_open( fname , 0);
    swat = ecl_file_iget_named_kw( ecl_file , "SWAT" , 0 );
    test_assert_true( util_double_approx_equal( ecl_kw_iget_float( swat , 0 ) , 1000 ));
  }
  test_work_area_free( work_area );
}
Ejemplo n.º 5
0
void test_grid(int nx, int ny, int nz) {
  ecl_kw_type * coord_kw = ecl_kw_alloc( COORD_KW , ECL_GRID_COORD_SIZE( nx , ny ) , ECL_FLOAT );
  ecl_kw_type * zcorn_kw = ecl_kw_alloc( ZCORN_KW , ECL_GRID_ZCORN_SIZE( nx , ny , nz) , ECL_FLOAT );
  int i,j,k;
  double a = 1.0;
  for (j= 0; j < ny; j++) {
    for (i = 0; i < nx; i++) {
      int offset = 6*(i + j*nx);
      ecl_kw_iset_float( coord_kw , offset    , a*i);
      ecl_kw_iset_float( coord_kw , offset + 1, a*j);
      ecl_kw_iset_float( coord_kw , offset + 2,  -1);

      ecl_kw_iset_float( coord_kw , offset + 3, a*i);
      ecl_kw_iset_float( coord_kw , offset + 4, a*j);
      ecl_kw_iset_float( coord_kw , offset + 5,  -1);

      for (k=0; k < nz; k++) {
        for (int c = 0; c < 4; c++) {
          int zi1 = ecl_grid_zcorn_index__( nx , ny , i , j , k , c);
          int zi2 = ecl_grid_zcorn_index__( nx , ny , i , j , k , c + 4);

          double z1 = k*a;
          double z2 = (k + 1) * a;

          ecl_kw_iset_float( zcorn_kw , zi1 , z1 );
          ecl_kw_iset_float( zcorn_kw , zi2 , z2 );
        }
      }
    }
  }


  {
    ecl_grid_type * grid = ecl_grid_alloc_GRDECL_kw( nx,ny,nz, zcorn_kw , coord_kw ,  NULL, NULL );
    test_assert_int_equal( ecl_grid_get_cell_twist3( grid , 0,0,0) , 0 );
    ecl_grid_free( grid );
  }


  {
    int zi1 = ecl_grid_zcorn_index__( nx , ny , 0 , 0 , 0 , 0);
    int zi2 = ecl_grid_zcorn_index__( nx , ny , 0 , 0 , 0 , 4);

    double z1 = 0;
    double z2 = -0.25;

    ecl_kw_iset_float( zcorn_kw , zi1 , z1 );
    ecl_kw_iset_float( zcorn_kw , zi2 , z2 );

    {
      ecl_grid_type * grid = ecl_grid_alloc_GRDECL_kw( nx,ny,nz, zcorn_kw , coord_kw ,  NULL, NULL );
      test_assert_int_equal( ecl_grid_get_cell_twist3( grid , 0,0,0) , 1 );
      ecl_grid_free( grid );
    }

    zi1 = ecl_grid_zcorn_index__( nx , ny , 0 , 0 , 0 , 3);
    zi2 = ecl_grid_zcorn_index__( nx , ny , 0 , 0 , 0 , 7);

    ecl_kw_iset_float( zcorn_kw , zi1 , z1 );
    ecl_kw_iset_float( zcorn_kw , zi2 , z2 );

    {
      ecl_grid_type * grid = ecl_grid_alloc_GRDECL_kw( nx,ny,nz, zcorn_kw , coord_kw ,  NULL, NULL );
      test_assert_int_equal( ecl_grid_get_cell_twist3( grid , 0,0,0) , 2 );
      ecl_grid_free( grid );
    }
  }


  ecl_kw_free( coord_kw );
  ecl_kw_free( zcorn_kw );
}
Ejemplo n.º 6
0
void ecl_rft_node_fwrite(const ecl_rft_node_type * rft_node, fortio_type * fortio, ert_ecl_unit_enum unit_set){
  ecl_rft_enum type = ecl_rft_node_get_type(rft_node);
  if (type != RFT)
    util_abort("%s: sorry - only writing of simple RFT is currently implemented",__func__);
  
  {
    ecl_kw_type * time = ecl_kw_alloc(TIME_KW, 1, ECL_FLOAT_TYPE);
    ecl_kw_iset_float(time, 0, ecl_rft_node_get_days(rft_node));
    ecl_kw_fwrite(time, fortio);
    ecl_kw_free(time);
  }

  {
    ecl_kw_type * datevalue = ecl_kw_alloc(DATE_KW, 3, ECL_INT_TYPE);
    time_t date = ecl_rft_node_get_date(rft_node);
    int day;
    int month;
    int year;
    ecl_util_set_date_values(date , &day , &month , &year);
    ecl_kw_iset_int(datevalue, 0, day);
    ecl_kw_iset_int(datevalue, 1, month);
    ecl_kw_iset_int(datevalue, 2, year);
    ecl_kw_fwrite(datevalue, fortio);
    ecl_kw_free(datevalue);
  }

  {
    ecl_kw_type * welletc = ecl_kw_alloc(WELLETC_KW, 16, ECL_CHAR_TYPE);
    ecl_rft_enum type = ecl_rft_node_get_type(rft_node);

    ecl_kw_iset_string8(welletc, 1, ecl_rft_node_get_well_name(rft_node));

    if(type == PLT) {
      ecl_kw_iset_string8(welletc, 5, "P");
    }else if(type == RFT){
      ecl_kw_iset_string8(welletc, 5, "R");
    }else if(type == SEGMENT){
      ecl_kw_iset_string8(welletc, 5, "S");
    }
    ecl_rft_node_fill_welletc(welletc, unit_set);
    ecl_kw_fwrite(welletc, fortio);
    ecl_kw_free(welletc);
  }

  {
    int size_cells = ecl_rft_node_get_size(rft_node);
    ecl_kw_type * conipos = ecl_kw_alloc(CONIPOS_KW, size_cells, ECL_INT_TYPE);
    ecl_kw_type * conjpos = ecl_kw_alloc(CONJPOS_KW, size_cells, ECL_INT_TYPE);
    ecl_kw_type * conkpos = ecl_kw_alloc(CONKPOS_KW, size_cells, ECL_INT_TYPE);
    ecl_kw_type * hostgrid = ecl_kw_alloc(HOSTGRID_KW, size_cells, ECL_CHAR_TYPE);
    ecl_kw_type * depth = ecl_kw_alloc(DEPTH_KW, size_cells, ECL_FLOAT_TYPE);
    ecl_kw_type * pressure = ecl_kw_alloc(PRESSURE_KW, size_cells, ECL_FLOAT_TYPE);
    ecl_kw_type * swat = ecl_kw_alloc(SWAT_KW, size_cells, ECL_FLOAT_TYPE);
    ecl_kw_type * sgas = ecl_kw_alloc(SGAS_KW, size_cells, ECL_FLOAT_TYPE);
    int i;

    for(i =0;i<size_cells;i++){
      const ecl_rft_cell_type * cell = vector_iget_const( rft_node->cells , i);
      ecl_kw_iset_int(conipos, i, ecl_rft_cell_get_i(cell)+1);
      ecl_kw_iset_int(conjpos, i, ecl_rft_cell_get_j(cell)+1);
      ecl_kw_iset_int(conkpos, i, ecl_rft_cell_get_k(cell)+1);
      ecl_kw_iset_float(depth, i, ecl_rft_cell_get_depth(cell));
      ecl_kw_iset_float(pressure, i, ecl_rft_cell_get_pressure(cell));
      ecl_kw_iset_float(swat, i, ecl_rft_cell_get_swat(cell));
      ecl_kw_iset_float(sgas, i, ecl_rft_cell_get_sgas(cell));
    }
    ecl_kw_fwrite(conipos, fortio);
    ecl_kw_fwrite(conjpos, fortio);
    ecl_kw_fwrite(conkpos, fortio);
    ecl_kw_fwrite(hostgrid, fortio);
    ecl_kw_fwrite(depth, fortio);
    ecl_kw_fwrite(pressure, fortio);
    ecl_kw_fwrite(swat, fortio);
    ecl_kw_fwrite(sgas, fortio);

    ecl_kw_free(conipos);
    ecl_kw_free(conjpos);
    ecl_kw_free(conkpos);
    ecl_kw_free(hostgrid);
    ecl_kw_free(depth);
    ecl_kw_free(pressure);
    ecl_kw_free(swat);
    ecl_kw_free(sgas);
  }
}