Beispiel #1
0
    void restoreRSandRV(const ecl_file_type* file_type,
                        SimulationConfigConstPtr sim_config,
                        int numcells,
                        BlackoilState* blackoil_state) {

        if (sim_config->hasDISGAS()) {
            const char* RS = "RS";
            if (ecl_file_has_kw(file_type, RS)) {
                ecl_kw_type* rs_kw = ecl_file_iget_named_kw(file_type, RS, 0);
                float* rs_data = ecl_kw_get_float_ptr(rs_kw);
                std::vector<double> rs_datavec(&rs_data[0], &rs_data[numcells]);
                blackoil_state->gasoilratio().clear();
                blackoil_state->gasoilratio().insert(blackoil_state->gasoilratio().begin(), rs_datavec.begin(), rs_datavec.end());
            } else {
                throw std::runtime_error("Restart file is missing RS data!\n");
            }
        }

        if (sim_config->hasVAPOIL()) {
            const char* RV = "RV";
            if (ecl_file_has_kw(file_type, RV)) {
                ecl_kw_type* rv_kw = ecl_file_iget_named_kw(file_type, RV, 0);
                float* rv_data = ecl_kw_get_float_ptr(rv_kw);
                std::vector<double> rv_datavec(&rv_data[0], &rv_data[numcells]);
                blackoil_state->rv().clear();
                blackoil_state->rv().insert(blackoil_state->rv().begin(), rv_datavec.begin(), rv_datavec.end());
            } else {
                throw std::runtime_error("Restart file is missing RV data!\n");
            }
        }
    }
Beispiel #2
0
    void restoreSaturation(const ecl_file_type* file_type,
                           const PhaseUsage& phaseUsage,
                           int numcells,
                           SimulatorState& simulator_state) {

        float* sgas_data = NULL;
        float* swat_data = NULL;

        if (phaseUsage.phase_used[BlackoilPhases::Aqua]) {
            const char* swat = "SWAT";
            if (ecl_file_has_kw(file_type, swat)) {
                ecl_kw_type* swat_kw = ecl_file_iget_named_kw(file_type , swat, 0);
                swat_data = ecl_kw_get_float_ptr(swat_kw);
                std::vector<double> swat_datavec(&swat_data[0], &swat_data[numcells]);
                EclipseIOUtil::addToStripedData(swat_datavec, simulator_state.saturation(), phaseUsage.phase_pos[BlackoilPhases::Aqua], phaseUsage.num_phases);
            } else {
                std::string error_str = "Restart file is missing SWAT data!\n";
                throw std::runtime_error(error_str);
            }
        }

        if (phaseUsage.phase_used[BlackoilPhases::Vapour]) {
            const char* sgas = "SGAS";
            if (ecl_file_has_kw(file_type, sgas)) {
                ecl_kw_type* sgas_kw = ecl_file_iget_named_kw(file_type , sgas, 0);
                sgas_data = ecl_kw_get_float_ptr(sgas_kw);
                std::vector<double> sgas_datavec(&sgas_data[0], &sgas_data[numcells]);
                EclipseIOUtil::addToStripedData(sgas_datavec, simulator_state.saturation(), phaseUsage.phase_pos[BlackoilPhases::Vapour], phaseUsage.num_phases);
            } else {
                std::string error_str = "Restart file is missing SGAS data!\n";
                throw std::runtime_error(error_str);
            }
        }
    }
Beispiel #3
0
int main(int argc , char ** argv) {
  const char * path_case = argv[1];
  char * grid_file = ecl_util_alloc_filename( NULL , path_case , ECL_EGRID_FILE , false , 0 );
  char * init_file = ecl_util_alloc_filename( NULL , path_case , ECL_INIT_FILE , false , 0 );

  ecl_file_type * init = ecl_file_open( init_file , 0 );
  ecl_grid_type * grid = ecl_grid_alloc( grid_file );
  const ecl_kw_type * poro_kw = ecl_file_iget_named_kw( init , "PORO" , 0 );
  const ecl_kw_type * porv_kw = ecl_file_iget_named_kw( init , "PORV" , 0 );
  ecl_kw_type * multpv = NULL;
  ecl_kw_type * NTG = NULL;
  bool error_found = false;

  double total_volume = 0;
  double total_diff = 0;
  int error_count = 0;
  int iactive;

  if (ecl_file_has_kw( init , "NTG"))
    NTG = ecl_file_iget_named_kw( init , "NTG" , 0);

  if (ecl_file_has_kw( init , "MULTPV"))
    multpv = ecl_file_iget_named_kw( init , "MULTPV" , 0);

  for (iactive = 0; iactive < ecl_grid_get_nactive( grid ); ++iactive) {
    int iglobal = ecl_grid_get_global_index1A( grid , iactive );
    double grid_volume = ecl_grid_get_cell_volume1( grid , iglobal );
    double eclipse_volume = ecl_kw_iget_float( porv_kw , iglobal ) / ecl_kw_iget_float( poro_kw , iactive );

    if (NTG)
      eclipse_volume /= ecl_kw_iget_float( NTG , iactive );

    if (multpv)
      eclipse_volume *= ecl_kw_iget_float( multpv , iactive);

    total_volume += grid_volume;
    total_diff += fabs( eclipse_volume - grid_volume );
    if (!util_double_approx_equal__( grid_volume , eclipse_volume , 2.5e-3, 0.00)) {
      double diff = 100 * (grid_volume - eclipse_volume) / eclipse_volume;
      printf("Error in cell: %d V1: %g    V2: %g   diff:%g %% \n", iglobal , grid_volume , eclipse_volume , diff);
      error_count++;
      error_found = true;
    }
  }
  printf("Total volume difference: %g %% \n", 100 * total_diff / total_volume );



  ecl_grid_free( grid );
  ecl_file_close( init );
  free( grid_file );
  free( init_file );
  if (error_found) {
    printf("Error_count: %d / %d \n",error_count , ecl_grid_get_nactive( grid ));
    exit(1);
  }  else
    exit(0);
}
Beispiel #4
0
bool well_state_add_MSW( well_state_type * well_state ,
                         const ecl_file_type * rst_file ,
                         int well_nr,
                         bool load_segment_information) {

  if (ecl_file_has_kw( rst_file , ISEG_KW)) {
    ecl_rsthead_type  * rst_head  = ecl_rsthead_alloc( rst_file );
    const ecl_kw_type * iwel_kw = ecl_file_iget_named_kw( rst_file , IWEL_KW , 0);
    const ecl_kw_type * iseg_kw = ecl_file_iget_named_kw( rst_file , ISEG_KW , 0);
    well_rseg_loader_type * rseg_loader = NULL;
    
    int segment_count;

    if (ecl_file_has_kw( rst_file , RSEG_KW )) {
      if (load_segment_information)
        rseg_loader = well_rseg_loader_alloc(rst_file);
        
      segment_count = well_segment_collection_load_from_kw( well_state->segments ,
                                                            well_nr ,
                                                            iwel_kw ,
                                                            iseg_kw ,
                                                            rseg_loader ,
                                                            rst_head,
                                                            load_segment_information , 
                                                            &well_state->is_MSW_well);
      
      
      if (segment_count > 0) {
        hash_iter_type * grid_iter = hash_iter_alloc( well_state->connections );
        while (!hash_iter_is_complete( grid_iter )) {
          const char * grid_name = hash_iter_get_next_key( grid_iter );
          const well_conn_collection_type * connections = hash_get( well_state->connections , grid_name );
          well_segment_collection_add_connections( well_state->segments , grid_name , connections );
        }
        hash_iter_free( grid_iter );

        well_segment_collection_link( well_state->segments );
        well_segment_collection_add_branches( well_state->segments , well_state->branches );
      }
      ecl_rsthead_free( rst_head );
      
      if (rseg_loader != NULL) {
        well_rseg_loader_free(rseg_loader);
      }
      
      return true;
    }
  } else
    return false;
}
Beispiel #5
0
static void ecl_rft_node_init_PLT_cells( ecl_rft_node_type * rft_node , const ecl_file_type * rft) {
  /* For PLT there is quite a lot of extra information which is not yet internalized. */
  const ecl_kw_type * conipos     = ecl_file_iget_named_kw( rft , CONIPOS_KW  , 0);
  const ecl_kw_type * conjpos     = ecl_file_iget_named_kw( rft , CONJPOS_KW  , 0);
  const ecl_kw_type * conkpos     = ecl_file_iget_named_kw( rft , CONKPOS_KW  , 0);  

  const int   * i      = ecl_kw_get_int_ptr( conipos );
  const int   * j      = ecl_kw_get_int_ptr( conjpos );
  const int   * k      = ecl_kw_get_int_ptr( conkpos );
  
  const float * WR               = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONWRAT_KW , 0));
  const float * GR               = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONGRAT_KW , 0)); 
  const float * OR               = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONORAT_KW , 0)); 
  const float * P                = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONPRES_KW , 0));
  const float * depth            = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONDEPTH_KW , 0));
  const float * flowrate         = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONVTUB_KW , 0));
  const float * oil_flowrate     = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONOTUB_KW , 0));
  const float * gas_flowrate     = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONGTUB_KW , 0));
  const float * water_flowrate   = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONWTUB_KW , 0));
  const float * connection_start = NULL;
  const float * connection_end   = NULL; 

  /* The keywords CONLENST_KW and CONLENEN_KW are ONLY present if we are dealing with a MSW well. */
  if (ecl_file_has_kw( rft , CONLENST_KW))
    connection_start = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONLENST_KW , 0));
  
  if (ecl_file_has_kw( rft , CONLENEN_KW))
    connection_end = ecl_kw_get_float_ptr( ecl_file_iget_named_kw( rft , CONLENEN_KW , 0));

  {
    int c;
    for ( c = 0; c < ecl_kw_get_size( conipos ); c++) {
      ecl_rft_cell_type * cell;
      double cs = 0;
      double ce = 0;
      
      if (connection_start)
        cs = connection_start[c];
      
      if (connection_end)
        ce = connection_end[c];

      /* The connection coordinates are shifted -= 1; i.e. all internal usage is offset 0. */
      cell = ecl_rft_cell_alloc_PLT( i[c] -1 , j[c] -1 , k[c] -1 , 
                                     depth[c] , P[c] , OR[c] , GR[c] , WR[c] , cs , ce,  flowrate[c] , oil_flowrate[c] , gas_flowrate[c] , water_flowrate[c]);
      ecl_rft_node_append_cell( rft_node , cell );
    }
  }
}
Beispiel #6
0
static void well_state_add_connections__( well_state_type * well_state ,
                                          const ecl_file_type * rst_file ,
                                          const char * grid_name ,
                                          int grid_nr,
                                          int well_nr ) {

  ecl_rsthead_type  * header   = ecl_rsthead_alloc( rst_file );
  const ecl_kw_type * icon_kw  = ecl_file_iget_named_kw( rst_file , ICON_KW   , 0);
  const ecl_kw_type * iwel_kw  = ecl_file_iget_named_kw( rst_file , IWEL_KW   , 0);

  
  well_state_add_wellhead( well_state , header , iwel_kw , well_nr , grid_name , grid_nr );

  if (!well_state_has_grid_connections( well_state , grid_name ))
    hash_insert_hash_owned_ref( well_state->connections , grid_name, well_conn_collection_alloc( ) , well_conn_collection_free__ );

  {
    ecl_kw_type * scon_kw = NULL;
    well_conn_collection_type * wellcc = hash_get( well_state->connections , grid_name );
    if (ecl_file_has_kw( rst_file , SCON_KW))
      scon_kw = ecl_file_iget_named_kw( rst_file , SCON_KW , 0);
    
    well_conn_collection_load_from_kw( wellcc , iwel_kw , icon_kw , scon_kw , well_nr , header );
  }
  ecl_rsthead_free( header );
}
Beispiel #7
0
static int well_state_get_lgr_well_nr( const well_state_type * well_state , const ecl_file_type * ecl_file) {
  int well_nr = -1;

  if (ecl_file_has_kw( ecl_file , ZWEL_KW)) {
    ecl_rsthead_type  * header  = ecl_rsthead_alloc( ecl_file );
    const ecl_kw_type * zwel_kw = ecl_file_iget_named_kw( ecl_file , ZWEL_KW  , 0 );
    int num_wells               = header->nwells;
    well_nr = 0;
    while (true) {
      bool found = false;
      {
        char * lgr_well_name = util_alloc_strip_copy( ecl_kw_iget_ptr( zwel_kw , well_nr * header->nzwelz) );

        if ( strcmp( well_state->name , lgr_well_name) == 0)
          found = true;
        else
          well_nr++;

        free( lgr_well_name );
      }

      if (found)
        break;
      else if (well_nr == num_wells) {
        // The well is not in this LGR at all.
        well_nr = -1;
        break;
      }

    }
    ecl_rsthead_free( header );
  }
  return well_nr;
}
Beispiel #8
0
    void restoreTemperatureData(const ecl_file_type* file,
                              EclipseStateConstPtr eclipse_state,
                              int numcells,
                              SimulatorState& simulator_state) {
        const char* temperature = "TEMP";

        if (ecl_file_has_kw(file , temperature)) {
            ecl_kw_type* temperature_kw = ecl_file_iget_named_kw(file, temperature, 0);

            if (ecl_kw_get_size(temperature_kw) != numcells) {
                throw std::runtime_error("Read of restart file: Could not restore temperature data, length of data from file not equal number of cells");
            }

            float* temperature_data = ecl_kw_get_float_ptr(temperature_kw);

            // factor and offset from the temperature values given in the deck to Kelvin
            double scaling = eclipse_state->getDeckUnitSystem().parse("Temperature")->getSIScaling();
            double offset  = eclipse_state->getDeckUnitSystem().parse("Temperature")->getSIOffset();

            for (size_t index = 0; index < simulator_state.temperature().size(); ++index) {
                simulator_state.temperature()[index] = unit::convert::from((double)temperature_data[index] - offset, scaling);
            }
          } else {
              throw std::runtime_error("Read of restart file: File does not contain TEMP data\n");
          }
    }
Beispiel #9
0
well_state_type * well_state_alloc_from_file( ecl_file_type * ecl_file , const ecl_grid_type * grid , int report_nr ,  int global_well_nr ,bool load_segment_information) {
  if (ecl_file_has_kw( ecl_file , IWEL_KW)) {
    well_state_type   * well_state = NULL;
    ecl_rsthead_type  * global_header  = ecl_rsthead_alloc( ecl_file );
    const ecl_kw_type * global_iwel_kw = ecl_file_iget_named_kw( ecl_file , IWEL_KW   , 0);
    const ecl_kw_type * global_zwel_kw = ecl_file_iget_named_kw( ecl_file , ZWEL_KW   , 0);

    const int iwel_offset = global_header->niwelz * global_well_nr;
    {
      char * name;
      bool open;
      well_type_enum type = UNDOCUMENTED_ZERO;
      {
        int int_state = ecl_kw_iget_int( global_iwel_kw , iwel_offset + IWEL_STATUS_ITEM );
        if (int_state > 0)
          open = true;
        else
          open = false;
      }

      {
        int int_type = ecl_kw_iget_int( global_iwel_kw , iwel_offset + IWEL_TYPE_ITEM);
        type = well_state_translate_ecl_type_int( int_type );
      }

      {
        const int zwel_offset         = global_header->nzwelz * global_well_nr;
        name = util_alloc_strip_copy(ecl_kw_iget_ptr( global_zwel_kw , zwel_offset ));  // Hardwired max 8 characters in Well Name
      }

      well_state = well_state_alloc(name , global_well_nr , open , type , report_nr , global_header->sim_time);
      free( name );

      well_state_add_connections( well_state , grid , ecl_file , global_well_nr);
      if (ecl_file_has_kw( ecl_file , ISEG_KW))
        well_state_add_MSW( well_state , ecl_file , global_well_nr , load_segment_information);
      

    }
    ecl_rsthead_free( global_header );
    return well_state;
  } else
    /* This seems a bit weird - have come over E300 restart files without the IWEL keyword. */
    return NULL;
}
Beispiel #10
0
int main(int argc , char ** argv) {
  test_install_SIGNALS();
  {
    const char * grid_file = argv[1];
    const char * rst_file_name = argv[2];

    ecl_grid_type * grid = ecl_grid_alloc( grid_file );
    ecl_file_type * rst_file = ecl_file_open( rst_file_name , 0);
    ecl_rsthead_type * header = ecl_rsthead_alloc( ecl_file_get_global_view( rst_file ) , ecl_util_filename_report_nr(rst_file_name) );
    const char * well_name = "WELL";
    int report_nr = 100;
    time_t valid_from = -1;
    bool open = false;
    well_type_enum type = ERT_GAS_INJECTOR;
    int global_well_nr = 0;
    bool load_segment_information = true;
    ecl_file_view_type * rst_view = ecl_file_get_global_view( rst_file );

    for (global_well_nr = 0; global_well_nr < header->nwells; global_well_nr++) {
      well_state_type * well_state = well_state_alloc(well_name , global_well_nr , open , type , report_nr , valid_from);
      test_assert_true( well_state_is_instance( well_state) );
      well_state_add_connections2( well_state , grid , rst_view , 0 );

      test_assert_true( well_state_has_grid_connections( well_state , ECL_GRID_GLOBAL_GRID ));
      test_assert_false( well_state_has_grid_connections( well_state , "???" ));
      test_assert_true( well_state_has_global_connections( well_state ));

      well_state_add_MSW2( well_state , rst_view , global_well_nr , load_segment_information );
      {
        const well_segment_collection_type * segments = well_state_get_segments( well_state );
        const well_branch_collection_type * branches = well_state_get_branches( well_state );

        if (well_state_is_MSW( well_state )) {
          test_assert_true( ecl_file_has_kw( rst_file , ISEG_KW ));
          test_assert_int_not_equal( well_segment_collection_get_size( segments ) , 0);
          test_assert_int_not_equal( well_branch_collection_get_size( branches ) , 0);
        } else {
          test_assert_int_equal( well_segment_collection_get_size( segments ) , 0);
          test_assert_int_equal( well_branch_collection_get_size( branches ) , 0);
          test_assert_false( well_state_is_MSW( well_state ));
        }
      }
      well_state_free( well_state );
    }
  }

  exit(0);
}
Beispiel #11
0
ecl_kw_type * ecl_nnc_export_get_tran_kw( const ecl_file_type * init_file , const char * kw , int lgr_nr ) {
  ecl_kw_type * tran_kw = NULL;
  if (lgr_nr == 0) {
    if (strcmp(kw , TRANNNC_KW) == 0)
      if(ecl_file_has_kw(init_file, kw)) {
        tran_kw = ecl_file_iget_named_kw(init_file, TRANNNC_KW, 0);
      }
  } else {
    if ((strcmp(kw , TRANNNC_KW) == 0) ||
        (strcmp(kw , TRANGL_KW) == 0)) {
      const int file_num_kw = ecl_file_get_size( init_file );
      int global_kw_index = 0;
      bool finished = false;
      bool correct_lgrheadi = false;
      int head_index = 0;
      int steps = 0;


      while(!finished){
        ecl_kw_type * ecl_kw = ecl_file_iget_kw( init_file , global_kw_index );
        const char *current_kw = ecl_kw_get_header(ecl_kw);
        if (strcmp( LGRHEADI_KW , current_kw) == 0) {
          if (ecl_kw_iget_int( ecl_kw , LGRHEADI_LGR_NR_INDEX) == lgr_nr) {
            correct_lgrheadi = true;
            head_index = global_kw_index;
          }else{
            correct_lgrheadi = false;
          }
        }
        if(correct_lgrheadi) {
          if (strcmp(kw, current_kw) == 0) {
            steps  = global_kw_index - head_index; /* This is to calculate who fare from lgrheadi we found the TRANGL/TRANNNC key word */
            if(steps == 3 || steps == 4 || steps == 6) { /* We only support a file format where TRANNNC is 3 steps and TRANGL is 4 or 6 steps from LGRHEADI */
              tran_kw = ecl_kw;
              finished = true;
              break;
            }
          }
        }
        global_kw_index++;
        if (global_kw_index == file_num_kw)
          finished = true;
      }
    }
  }
  return tran_kw;
}
static ecl_grav_survey_type * ecl_grav_survey_alloc_RPORV(ecl_grav_type * ecl_grav , 
                                                          const ecl_file_type * restart_file , 
                                                          const char * name ) {
  ecl_grav_survey_type * survey = ecl_grav_survey_alloc_empty( ecl_grav , name , GRAV_CALC_RPORV);
  if (ecl_file_has_kw( restart_file , RPORV_KW)) {
    ecl_kw_type * rporv_kw = ecl_file_iget_named_kw( restart_file , RPORV_KW , 0);
    int iactive;
    for (iactive = 0; iactive < ecl_kw_get_size( rporv_kw ); iactive++) 
      survey->porv[ iactive ] = ecl_kw_iget_as_double( rporv_kw , iactive );
  } else 
    util_abort("%s: restart file did not contain %s keyword??\n",__func__ , RPORV_KW);
  
  {
    const ecl_file_type * init_file = ecl_grav->init_file;
    ecl_grav_survey_assert_RPORV( survey , init_file );
    ecl_grav_survey_add_phases( ecl_grav , survey ,  restart_file , GRAV_CALC_RPORV);
  }
  return survey;
}
//--------------------------------------------------------------------------------------------------
/// 
//--------------------------------------------------------------------------------------------------
ecl_kw_type* RifEclipseOutputFileTools::createActnumFromPorv(ecl_file_type* ecl_file)
{
    std::string porv_kw("PORV");

    if (ecl_file_has_kw(ecl_file, porv_kw.data()))
    {
        ecl_file_view_type* fileView = ecl_file_get_global_view(ecl_file);

        int keywordCount = ecl_file_get_num_named_kw(ecl_file, porv_kw.data());
        for (int index = 0; index < keywordCount; index++)
        {
            ecl_kw_type* fileKeyword = ecl_file_view_iget_named_kw(fileView, porv_kw.data(), index);
            if (fileKeyword)
            {
                float porvLimit = 0.0f;

                return ecl_kw_alloc_actnum(fileKeyword, porvLimit);
            }
        }
    }

    return nullptr;
}
Beispiel #14
0
    void restorePressureData(const ecl_file_type* file,
                             EclipseStateConstPtr eclipse_state,
                             int numcells,
                             SimulatorState& simulator_state) {
        const char* pressure = "PRESSURE";

        if (ecl_file_has_kw(file , pressure)) {

            ecl_kw_type* pressure_kw = ecl_file_iget_named_kw(file, pressure, 0);

            if (ecl_kw_get_size(pressure_kw) != numcells) {
                throw std::runtime_error("Read of restart file: Could not restore pressure data, length of data from file not equal number of cells");
            }

            float* pressure_data = ecl_kw_get_float_ptr(pressure_kw);
            const double deck_pressure_unit = (eclipse_state->getDeckUnitSystem().getType() == UnitSystem::UNIT_TYPE_METRIC) ? Opm::unit::barsa : Opm::unit::psia;
            for (size_t index = 0; index < simulator_state.pressure().size(); ++index) {
                simulator_state.pressure()[index] = unit::convert::from((double)pressure_data[index], deck_pressure_unit);
            }
        } else {
            throw std::runtime_error("Read of restart file: File does not contain PRESSURE data\n");
        }
    }
Beispiel #15
0
ecl_rft_node_type * ecl_rft_node_alloc(const ecl_file_type * rft) {
  ecl_kw_type       * welletc   = ecl_file_iget_named_kw(rft , WELLETC_KW , 0);
  ecl_rft_node_type * rft_node  = ecl_rft_node_alloc_empty(ecl_kw_iget_ptr(welletc , WELLETC_TYPE_INDEX));
  
  if (rft_node != NULL) {
    ecl_kw_type * date_kw = ecl_file_iget_named_kw( rft , DATE_KW    , 0);
    rft_node->well_name = util_alloc_strip_copy( ecl_kw_iget_ptr(welletc , WELLETC_NAME_INDEX));
    
    /* Time information. */
    {
      int * time = ecl_kw_get_int_ptr( date_kw );
      rft_node->recording_date = ecl_util_make_date( time[DATE_DAY_INDEX] , time[DATE_MONTH_INDEX] , time[DATE_YEAR_INDEX] );
    }
    rft_node->days = ecl_kw_iget_float( ecl_file_iget_named_kw( rft , TIME_KW , 0 ) , 0);
    if (ecl_file_has_kw( rft , CONLENST_KW))
      rft_node->MSW = true;
    else
      rft_node->MSW = false;

    ecl_rft_node_init_cells( rft_node , rft );
  }
  return rft_node;
}
Beispiel #16
0
static double gravity_response(const ecl_grid_type * ecl_grid      , 
                               const ecl_file_type * init_file     , 
                               const ecl_file_type * restart_file1 , 
                               const ecl_file_type * restart_file2 ,
                               const grav_station_type * grav_station , 
                               int model_phases, 
                               int file_phases) {
  
  ecl_kw_type * rporv1_kw   = NULL;  
  ecl_kw_type * rporv2_kw   = NULL;
  ecl_kw_type * oil_den1_kw = NULL;  
  ecl_kw_type * oil_den2_kw = NULL;
  ecl_kw_type * gas_den1_kw = NULL;
  ecl_kw_type * gas_den2_kw = NULL;
  ecl_kw_type * wat_den1_kw = NULL;
  ecl_kw_type * wat_den2_kw = NULL;
  ecl_kw_type * sgas1_kw    = NULL;
  ecl_kw_type * sgas2_kw    = NULL;
  ecl_kw_type * swat1_kw    = NULL;
  ecl_kw_type * swat2_kw    = NULL;
  ecl_kw_type * aquifern_kw = NULL ;
  double local_deltag = 0;

  /* Extracting the pore volumes */
  rporv1_kw = ecl_file_iget_named_kw( restart_file1 , "RPORV" , 0);      
  rporv2_kw = ecl_file_iget_named_kw( restart_file2 , "RPORV" , 0);      
  
  
  /** Extracting the densities */
  {
    // OIL_DEN
    if( has_phase(model_phases , OIL) ) {
      if (simulator == ECLIPSE100) {
        oil_den1_kw  = ecl_file_iget_named_kw(restart_file1, "OIL_DEN", 0);
        oil_den2_kw  = ecl_file_iget_named_kw(restart_file2, "OIL_DEN", 0);
      } else { // ECLIPSE300
        oil_den1_kw  = ecl_file_iget_named_kw(restart_file1, "DENO", 0);
        oil_den2_kw  = ecl_file_iget_named_kw(restart_file2, "DENO", 0);
      } ;
    }
    
    // GAS_DEN
    if( has_phase( model_phases , GAS) ) {
      if (simulator == ECLIPSE100) {
        gas_den1_kw  = ecl_file_iget_named_kw(restart_file1, "GAS_DEN", 0);
        gas_den2_kw  = ecl_file_iget_named_kw(restart_file2, "GAS_DEN", 0);
      } else { // ECLIPSE300
        gas_den1_kw  = ecl_file_iget_named_kw(restart_file1, "DENG", 0);
        gas_den2_kw  = ecl_file_iget_named_kw(restart_file2, "DENG", 0);
      } ;
    }
    
    // WAT_DEN
    if( has_phase( model_phases , WATER) ) {
      if (simulator == ECLIPSE100) {
        wat_den1_kw  = ecl_file_iget_named_kw(restart_file1, "WAT_DEN", 0);
        wat_den2_kw  = ecl_file_iget_named_kw(restart_file2, "WAT_DEN", 0);
      } else { // ECLIPSE300
        wat_den1_kw  = ecl_file_iget_named_kw(restart_file1, "DENW", 0);
        wat_den2_kw  = ecl_file_iget_named_kw(restart_file2, "DENW", 0);
      } ;
    }
  }
  
  
  /* Extracting the saturations */
  {
    // SGAS
    if( has_phase( file_phases , GAS )) {
      sgas1_kw     = ecl_file_iget_named_kw(restart_file1, "SGAS", 0);
      sgas2_kw     = ecl_file_iget_named_kw(restart_file2, "SGAS", 0);
    } 
    
    // SWAT
    if( has_phase( file_phases , WATER )) {
      swat1_kw     = ecl_file_iget_named_kw(restart_file1, "SWAT", 0);
      swat2_kw     = ecl_file_iget_named_kw(restart_file2, "SWAT", 0);
    } 
  }
  
  
  /* The numerical aquifer information */
  if( ecl_file_has_kw( init_file , "AQUIFERN")) 
    aquifern_kw     = ecl_file_iget_named_kw(init_file, "AQUIFERN", 0);
  {
    int     nactive  = ecl_grid_get_active_size( ecl_grid );
    float * zero     = util_calloc( nactive , sizeof * zero     );    /* Fake vector of zeros used for densities / sturations when you do not have data. */
    int   * int_zero = util_calloc( nactive , sizeof * int_zero );    /* Fake vector of zeros used for AQUIFER when the init file does not supply data. */
    /* 
       Observe that the fake vectors are only a coding simplification,
       they should not be really used.
    */

    {
      int i;
      for (i=0; i < nactive; i++) {
        zero[i]     = 0;
        int_zero[i] = 0;
      }
    }
    {
      const float * sgas1_v   = safe_get_float_ptr( sgas1_kw    , NULL );
      const float * swat1_v   = safe_get_float_ptr( swat1_kw    , NULL );
      const float * oil_den1  = safe_get_float_ptr( oil_den1_kw , zero );
      const float * gas_den1  = safe_get_float_ptr( gas_den1_kw , zero );
      const float * wat_den1  = safe_get_float_ptr( wat_den1_kw , zero );
      
      const float * sgas2_v   = safe_get_float_ptr( sgas2_kw    , NULL );
      const float * swat2_v   = safe_get_float_ptr( swat2_kw    , NULL );
      const float * oil_den2  = safe_get_float_ptr( oil_den2_kw , zero );
      const float * gas_den2  = safe_get_float_ptr( gas_den2_kw , zero );
      const float * wat_den2  = safe_get_float_ptr( wat_den2_kw , zero );
      
      const float * rporv1    = ecl_kw_get_float_ptr(rporv1_kw);
      const float * rporv2    = ecl_kw_get_float_ptr(rporv2_kw);
      double utm_x = grav_station->utm_x;
      double utm_y = grav_station->utm_y;
      double tvd   = grav_station->depth;
      
      int   * aquifern;
      int global_index;
          
      if (aquifern_kw != NULL)
        aquifern = ecl_kw_get_int_ptr( aquifern_kw );
      else
        aquifern = int_zero;

      for (global_index=0;global_index < ecl_grid_get_global_size( ecl_grid ); global_index++){
        const int act_index = ecl_grid_get_active_index1( ecl_grid , global_index );
        if (act_index >= 0) {

          // Not numerical aquifer 
          if(aquifern[act_index] >= 0){ 
            float swat1 = swat1_v[act_index];
            float swat2 = swat2_v[act_index];
            float sgas1 = 0;
            float sgas2 = 0;
            float soil1 = 0;
            float soil2 = 0;

            truncate_saturation( &swat1 );
            truncate_saturation( &swat2 );
            
            if (has_phase( model_phases , GAS)) {
              if (has_phase( file_phases , GAS )) {
                sgas1 = sgas1_v[act_index];
                sgas2 = sgas2_v[act_index];
                truncate_saturation( &sgas1 );
                truncate_saturation( &sgas2 );
              } else {
                sgas1 = 1 - swat1;
                sgas2 = 1 - swat2;
              }
            }
            
            if (has_phase( model_phases , OIL )) {
              soil1 =  1 - sgas1  - swat1;
              soil2 =  1 - sgas2  - swat2;
              truncate_saturation( &soil1 );
              truncate_saturation( &soil2 );
            }
            
                        
            /* 
               We have found all the info we need for one cell.
            */
            
            {
              double  mas1 , mas2;
              double  xpos , ypos , zpos;
              
              mas1 = rporv1[act_index]*(soil1 * oil_den1[act_index] + sgas1 * gas_den1[act_index] + swat1 * wat_den1[act_index] );
              mas2 = rporv2[act_index]*(soil2 * oil_den2[act_index] + sgas2 * gas_den2[act_index] + swat2 * wat_den2[act_index] );
              
              ecl_grid_get_xyz1(ecl_grid , global_index , &xpos , &ypos , &zpos);
              {
                double dist_x   = xpos - utm_x;
                double dist_y   = ypos - utm_y;
                double dist_d   = zpos - tvd;
                double dist_sq  = dist_x*dist_x + dist_y*dist_y + dist_d*dist_d;
                
                if(dist_sq == 0){
                  exit(1);
                }
                local_deltag += 6.67428E-3*(mas2 - mas1)*dist_d/pow(dist_sq, 1.5); // Gravity in units of \mu Gal = 10^{-8} m/s^2
              }
              
            }
          }
        }
      }
    }
    free( zero );
    free( int_zero );
  }
  return local_deltag;
}
Beispiel #17
0
static int gravity_check_input( const ecl_grid_type * ecl_grid , 
                                const ecl_file_type * init_file , 
                                const ecl_file_type * restart_file1, 
                                const ecl_file_type * restart_file2,
                                int   * __model_phases,
                                int   * __file_phases) {
  {
    int model_phases = 0;
    int file_phases  = 0;

    /* Check which phases are present in the model */
    if (ecl_file_has_kw(restart_file1 , "OIL_DEN")) {
      model_phases += OIL;  
      simulator = ECLIPSE100 ;
    } else if (ecl_file_has_kw(restart_file1 , "DENO")) {
      model_phases += OIL;  
      simulator = ECLIPSE300 ;
    } ;
      
    if (ecl_file_has_kw(restart_file1 , "WAT_DEN")) {
      model_phases += WATER;                         
      simulator = ECLIPSE100 ;
    } else if (ecl_file_has_kw(restart_file1 , "DENW")) {
      model_phases += WATER;                         
      simulator = ECLIPSE300 ;
    } ;
    
    if (ecl_file_has_kw(restart_file1 , "GAS_DEN")) {
      model_phases += GAS;
      simulator = ECLIPSE100 ;
    } else if (ecl_file_has_kw(restart_file1 , "DENG")) {
      model_phases += GAS;
      simulator = ECLIPSE300 ;
    } ;
    
    
    /* Check which phases are present in the restart files. We assume the restart file NEVER has SOIL information */
    if (ecl_file_has_kw(restart_file1 , "SWAT"))
      file_phases += WATER;
    if (ecl_file_has_kw(restart_file1 , "SGAS"))
      file_phases += GAS;
    
    
    /* Consiency check */
    {
      /**
         The following assumptions are made:
         
         1. All restart files should have water, i.e. the SWAT keyword. 
         2. All phases present in the restart file should also be present as densities, 
            in addition the model must contain one additional phase. 
         3. The restart files can never contain oil saturation.
         
      */
      if ( !has_phase( file_phases , WATER ) )
        util_exit("Could not locate SWAT keyword in restart files\n");
      
      if ( has_phase( file_phases , OIL ))
        util_exit("Can not handle restart files with SOIL keyword\n"); 
      
      if (! has_phase( model_phases , WATER ) )
        util_exit("Could not locate WAT_DEN keyword in restart files\n");      
      
      if ( has_phase( file_phases , GAS )) {
        /** Restart file has both water and gas - means we need all three densities. */
        if (! (has_phase( model_phases , GAS) && has_phase( model_phases , OIL)))
          util_exit("Could not find GAS_DEN and OIL_DEN keywords in restart files\n");
      } else {
        /* This is (water + oil) or (water + gas) system. We enforce one of the densities.*/
        if ( !has_phase( model_phases , GAS + OIL))
          util_exit("Could not find either GAS_DEN or OIL_DEN kewyords in restart files\n");
      }
    }
    *__model_phases = model_phases;
    *__file_phases  = file_phases;
  }
  
  /* Check that the restart files have RPORV information. This is ensured by giving the argument RPORV to the RPTRST keyword. */
  if ( !(ecl_file_has_kw( restart_file1 , "RPORV") && ecl_file_has_kw( restart_file2 , "RPORV")) )
    util_exit("Sorry: the restartfiles do  not contain RPORV\n");       


  /**
     Check that the rporv values are in the right ballpark.  For
     ECLIPSE version 2008.2 they are way f*****g off. Check PORV
     versus RPORV for ten 'random' locations in the grid.
  */
  {
    const ecl_kw_type * rporv1_kw     = ecl_file_iget_named_kw( restart_file1 , "RPORV" , 0);      
    const ecl_kw_type * rporv2_kw     = ecl_file_iget_named_kw( restart_file2 , "RPORV" , 0);      
    const ecl_kw_type * init_porv_kw  = ecl_file_iget_named_kw( init_file     , "PORV" , 0);

    int    active_index;
    int    active_delta;
    int    active_size;
    
    ecl_grid_get_dims( ecl_grid , NULL , NULL , NULL , &active_size );
    active_delta = active_size / 12;
    for (active_index = active_delta; active_index < active_size; active_index += active_delta) {
      int    global_index = ecl_grid_get_global_index1A( ecl_grid , active_index );
      double init_porv    = ecl_kw_iget_as_double( init_porv_kw , global_index );   /* NB - this uses global indexing. */
      double rporv1       = ecl_kw_iget_as_double( rporv1_kw ,  active_index );
      double rporv2       = ecl_kw_iget_as_double( rporv2_kw ,  active_index );
      double rporv12      = 0.5 * ( rporv1 + rporv2 );
      double fraction     = util_double_min( init_porv , rporv12 ) / util_double_max( init_porv , rporv12 );

      if (fraction  < 0.50) {
        fprintf(stderr,"-----------------------------------------------------------------\n");
        fprintf(stderr,"INIT PORV: %g \n",init_porv);
        fprintf(stderr,"RPORV1   : %g \n",rporv1);
        fprintf(stderr,"RPORV2   : %g \n",rporv2);
        fprintf(stderr,"Hmmm - the RPORV values extracted from the restart file seem to be \n");
        fprintf(stderr,"veeery different from the initial rporv value. This might indicated\n");
        fprintf(stderr,"an ECLIPSE bug. Version 2007.2 is known to be ok in this respect, \n");
        fprintf(stderr,"whereas version 2008.2 is known to have a bug. \n");
        fprintf(stderr,"-----------------------------------------------------------------\n");
        exit(1);
      }
    }
  }

  return 0;
}
static ecl_grav_phase_type * ecl_grav_phase_alloc( ecl_grav_type * ecl_grav , 
                                                   ecl_grav_survey_type * survey , 
                                                   ecl_phase_enum phase , 
                                                   const ecl_file_type * restart_file, 
                                                   grav_calc_type calc_type) {
  
  const ecl_file_type * init_file        = ecl_grav->init_file;
  const ecl_grid_cache_type * grid_cache = ecl_grav->grid_cache;
  const char * sat_kw_name               = ecl_util_get_phase_name( phase );
  {
    ecl_grav_phase_type * grav_phase = util_malloc( sizeof * grav_phase );
    const int size                   = ecl_grid_cache_get_size( grid_cache );
    
    UTIL_TYPE_ID_INIT( grav_phase , ECL_GRAV_PHASE_TYPE_ID );
    grav_phase->grid_cache   = grid_cache;
    grav_phase->aquifer_cell = ecl_grav->aquifer_cell;
    grav_phase->fluid_mass   = util_calloc( size , sizeof * grav_phase->fluid_mass );
    grav_phase->phase        = phase;
    grav_phase->work         = NULL;

    if (calc_type == GRAV_CALC_FIP) {
      ecl_kw_type * pvtnum_kw = ecl_file_iget_named_kw( init_file , PVTNUM_KW , 0 );
      double_vector_type * std_density = hash_get( ecl_grav->std_density , ecl_util_get_phase_name( phase ));
      ecl_kw_type * fip_kw;
        
      if ( phase == ECL_OIL_PHASE)
        fip_kw = ecl_file_iget_named_kw( restart_file , FIPOIL_KW , 0 );
      else if (phase == ECL_GAS_PHASE)
        fip_kw = ecl_file_iget_named_kw( restart_file , FIPGAS_KW , 0 );
      else
        fip_kw = ecl_file_iget_named_kw( restart_file , FIPWAT_KW , 0 );
      
      {
        int iactive;
        for (iactive=0; iactive < size; iactive++) {
          double fip    = ecl_kw_iget_as_double( fip_kw , iactive );
          int    pvtnum = ecl_kw_iget_int( pvtnum_kw , iactive );
          grav_phase->fluid_mass[ iactive ] = fip * double_vector_safe_iget( std_density , pvtnum );
        }
      }
    } else {
      ecl_version_enum      ecl_version = ecl_file_get_ecl_version( init_file );
      const char          * den_kw_name = get_den_kw( phase , ecl_version );
      const ecl_kw_type   * den_kw      = ecl_file_iget_named_kw( restart_file , den_kw_name , 0 );

      if (calc_type == GRAV_CALC_RFIP) {
        ecl_kw_type * rfip_kw;
        if ( phase == ECL_OIL_PHASE)
          rfip_kw = ecl_file_iget_named_kw( restart_file , RFIPOIL_KW , 0 );
        else if (phase == ECL_GAS_PHASE)
          rfip_kw = ecl_file_iget_named_kw( restart_file , RFIPGAS_KW , 0 );
        else
          rfip_kw = ecl_file_iget_named_kw( restart_file , RFIPWAT_KW , 0 );
        
        {
          int iactive;
          for (iactive=0; iactive < size; iactive++) {
            double rho   = ecl_kw_iget_as_double( den_kw  , iactive );
            double rfip  = ecl_kw_iget_as_double( rfip_kw , iactive );
            grav_phase->fluid_mass[ iactive ] = rho * rfip;
          }
        }
      } else {
        /* (calc_type == GRAV_CALC_RPORV) || (calc_type == GRAV_CALC_PORMOD) */
        ecl_kw_type * sat_kw;
        bool private_sat_kw = false;
        if (ecl_file_has_kw( restart_file , sat_kw_name )) 
          sat_kw = ecl_file_iget_named_kw( restart_file , sat_kw_name , 0 );
        else {
          /* We are targeting the residual phase, e.g. the OIL phase in a three phase system. */
          const ecl_kw_type * swat_kw = ecl_file_iget_named_kw( restart_file , "SWAT" , 0 );
          sat_kw = ecl_kw_alloc_copy( swat_kw );
          ecl_kw_scalar_set_float( sat_kw , 1.0 );
          ecl_kw_inplace_sub( sat_kw , swat_kw );  /* sat = 1 - SWAT */
          
          if (ecl_file_has_kw( restart_file , "SGAS" )) {
            const ecl_kw_type * sgas_kw = ecl_file_iget_named_kw( restart_file , "SGAS" , 0 );
            ecl_kw_inplace_sub( sat_kw , sgas_kw );  /* sat -= SGAS */
          }
          private_sat_kw = true;
        }
        
        {
          int iactive;
          for (iactive=0; iactive < size; iactive++) {
            double rho  = ecl_kw_iget_as_double( den_kw , iactive );
            double sat  = ecl_kw_iget_as_double( sat_kw , iactive );
            grav_phase->fluid_mass[ iactive ] = rho * sat * survey->porv[ iactive ];
          }
        }
        
        if (private_sat_kw)
          ecl_kw_free( sat_kw );
      }
    }
    
    return grav_phase;
  }
}