Skip to content
Snippets Groups Projects
Field.cpp 33.2 KiB
Newer Older
// $Id: Field.cpp,v 1.39 2008-06-12 09:31:36 geuzaine Exp $
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed
//
// Copyright (C) 1997-2008 C. Geuzaine, J.-F. Remacle
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed
//
// This program is free software; you can redistribute it and/or modify
// it under the terms of the GNU General Public License as published by
// the Free Software Foundation; either version 2 of the License, or
// (at your option) any later version.
//
// This program is distributed in the hope that it will be useful,
// but WITHOUT ANY WARRANTY; without even the implied warranty of
// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU General Public License
// along with this program; if not, write to the Free Software
// Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307
// USA.
// 
// Please report all bugs and problems to <gmsh@geuz.org>.

#include <list>
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed
#include <math.h>
#include <fstream>
#include <string>
#include <sstream>

#ifdef HAVE_MATH_EVAL
#include "matheval.h"
#endif
#ifdef HAVE_ANN
#include "ANN/ANN.h"
#endif

#include "Context.h"
#include "Field.h"
#include "GeoInterpolation.h"
#include "GModel.h"
#include "Message.h"
#include "OctreePost.h"
#include "PViewDataList.h"
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed

#define MAX_LC 1.e22

extern Context_T CTX;
class FieldOptionDouble:public FieldOption
{
public:
Jean-François Remacle's avatar
Jean-François Remacle committed
  double &val;
  FieldOptionType get_type()
  {
    return FIELD_OPTION_DOUBLE;
  };
FieldOptionDouble(double &_val, bool * _status = NULL):FieldOption(_status),
    val(_val) {
  };
Jean-François Remacle's avatar
Jean-François Remacle committed
  double numerical_value() const
  {
    return val;
Jean-François Remacle's avatar
Jean-François Remacle committed
  void numerical_value(double v)
  {
    modified();
    val = v;
Jean-François Remacle's avatar
Jean-François Remacle committed
  void get_text_representation(std::string & v_str)
  {
    std::ostringstream sstream;
    sstream.precision(16);
    sstream << val;
    v_str = sstream.str();
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  int &val;
  FieldOptionType get_type()
  {
    return FIELD_OPTION_INT;
  };
FieldOptionInt(int &_val, bool * _status = NULL):FieldOption(_status),
    val(_val) {
  };
Jean-François Remacle's avatar
Jean-François Remacle committed
  double numerical_value() const
  {
    return val;
Jean-François Remacle's avatar
Jean-François Remacle committed
  void numerical_value(double v)
  {
    modified();
Ruth Sabariego's avatar
Ruth Sabariego committed
    val = (int)v;
Jean-François Remacle's avatar
Jean-François Remacle committed
  void get_text_representation(std::string & v_str)
  {
    std::ostringstream sstream;
    sstream << val;
    v_str = sstream.str();
  }
class FieldOptionList:public FieldOption
{
public:
  std::list < int >&val;
Jean-François Remacle's avatar
Jean-François Remacle committed
  FieldOptionType get_type()
  {
    return FIELD_OPTION_LIST;
  };
FieldOptionList(std::list < int >&_val, bool * _status = NULL):FieldOption(_status),
    val(_val)
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
    modified();
    return val;
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
    return val;
  }
  void get_text_representation(std::string & v_str)
  {
    std::ostringstream sstream;
    sstream << "{";
    for(std::list < int >::iterator it = val.begin(); it != val.end(); it++) {
      if(it != val.begin())
        sstream << ", ";
      sstream << *it;
    }
    sstream << "}";
    v_str = sstream.str();
  }
class FieldOptionString:public FieldOption
{
public:
  std::string & val;
Jean-François Remacle's avatar
Jean-François Remacle committed
  virtual FieldOptionType get_type()
  {
    return FIELD_OPTION_STRING;
  };
FieldOptionString(std::string & _val, bool * _status = NULL):FieldOption(_status),
    val(_val) {
  };
  std::string & string() {
Jean-François Remacle's avatar
Jean-François Remacle committed
    modified();
    return val;
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
    return val;
  }
  void get_text_representation(std::string & v_str)
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
    std::ostringstream sstream;
    sstream << "\"" << val << "\"";
    v_str = sstream.str();
  }
class FieldOptionBool:public FieldOption
{
public:
  bool & val;
Jean-François Remacle's avatar
Jean-François Remacle committed
  FieldOptionType get_type()
  {
    return FIELD_OPTION_BOOL;
  };
FieldOptionBool(bool & _val, bool * _status = NULL):FieldOption(_status),
    val(_val) {
  };
Jean-François Remacle's avatar
Jean-François Remacle committed
  double numerical_value() const
  {
    return val;
Jean-François Remacle's avatar
Jean-François Remacle committed
  void numerical_value(double v)
  {
    modified();
    val = v;
Jean-François Remacle's avatar
Jean-François Remacle committed
  void get_text_representation(std::string & v_str)
  {
    std::ostringstream sstream;
    sstream << val;
    v_str = sstream.str();
  }
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed
void FieldManager::reset()
{
  for(std::map < int, Field * >::iterator it = begin(); it != end(); it++) {
    delete it->second;
  }
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed
Field *FieldManager::get(int id)
{
  iterator it = find(id);
Jean-François Remacle's avatar
Jean-François Remacle committed
  if(it == end()) {
    return NULL;
  }
  return it->second;
}
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed

Field *FieldManager::new_field(int id, std::string type_name)
Jean-François Remacle's avatar
Jean-François Remacle committed
  if(find(id) != end()) {
    Msg::Error("Field id %i is already defined.", id);
Jean-François Remacle's avatar
Jean-François Remacle committed
    return NULL;
  }
  if(map_type_name.find(type_name) == map_type_name.end()) {
    Msg::Error("Unknown field type \"%s\".", type_name.c_str());
Jean-François Remacle's avatar
Jean-François Remacle committed
    return NULL;
  }
  Field *f = (*map_type_name[type_name]) ();
  if(!f)
    return NULL;
  f->id = id;
  (*this)[id] = f;
Jean-François Remacle's avatar
Jean-François Remacle committed
int FieldManager::new_id()
{
  int i = 0;
  iterator it = begin();
  while(1) {
    i++;
    while(it != end() && it->first < i)
      it++;
    if(it == end() || it->first != i)
      break;
  }
  return std::max(i, 1);
Jean-François Remacle's avatar
Jean-François Remacle committed

int FieldManager::max_id()
{
  if(!empty())
    return rbegin()->first;
  else
    return 0;
Jean-François Remacle's avatar
Jean-François Remacle committed
void FieldManager::delete_field(int id)
{
  iterator it = find(id);
  if(it == end()) {
    Msg::Error("Cannot delete field id %i, it does not exist.", id);
Jean-François Remacle's avatar
Jean-François Remacle committed
    return;
  }
  delete it->second;
  erase(it);
Christophe Geuzaine's avatar
 
Christophe Geuzaine committed
// StructuredField
  double o[3], d[3];
  int n[3];
  double *data;
Jean-François Remacle's avatar
Jean-François Remacle committed
  bool error_status;
Jean-François Remacle's avatar
Jean-François Remacle committed
  std::string file_name;
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
    options["FileName"] = new FieldOptionString(file_name, &update_needed);
    text_format = false;
    options["TextFormat"] = new FieldOptionBool(text_format, &update_needed);
Jean-François Remacle's avatar
Jean-François Remacle committed
    data = NULL;
  }
  const char *get_name()
  {
    return "Structured";
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  }
  double operator() (double x, double y, double z)
  {
    if(update_needed) {
      error_status = false;
      try {
        std::ifstream input(file_name.c_str());
        if(!input.is_open())
          throw(1);
        input.
          exceptions(std::ifstream::eofbit | std::ifstream::failbit | std::
                     ifstream::badbit);
        if(!text_format) {
          input.read((char *)o, 3 * sizeof(double));
          input.read((char *)d, 3 * sizeof(double));
          input.read((char *)n, 3 * sizeof(int));
          int nt = n[0] * n[1] * n[2];
          if(data)
            delete[]data;
          data = new double[nt];
          input.read((char *)data, nt * sizeof(double));
        }
        else {
          input >> o[0] >> o[1] >> o[2] >> d[0] >> d[1] >> d[2] >> n[0] >>
            n[1] >> n[2];
          int nt = n[0] * n[1] * n[2];
          if(data)
            delete[]data;
          data = new double[nt];
          for(int i = 0; i < nt; i++)
            input >> data[i];
        }
Jean-François Remacle's avatar
Jean-François Remacle committed
        input.close();
      }
      catch(...) {
        error_status = true;
        Msg::Error("Field %i : error reading file %s", this->id,
Jean-François Remacle's avatar
Jean-François Remacle committed
            file_name.c_str());
      }
      update_needed = false;
    }
    if(error_status)
      return MAX_LC;
    //tri-linear
    int id[2][3];
    double xi[3];
    double xyz[3] = { x, y, z };
    for(int i = 0; i < 3; i++) {
      id[0][i] = (int)floor((xyz[i] - o[i]) / d[i]);
      id[1][i] = id[0][i] + 1;
      id[0][i] = std::max(std::min(id[0][i], n[i] - 1), 0);
      id[1][i] = std::max(std::min(id[1][i], n[i] - 1), 0);
      xi[i] = (xyz[i] - (o[i] + id[0][i] * d[i])) / d[i];
      xi[i] = std::max(std::min(xi[i], 1.), 0.);
    }
    double v = 0;
    for(int i = 0; i < 2; i++)
      for(int j = 0; j < 2; j++)
        for(int k = 0; k < 2; k++) {
          v += data[id[i][0] * n[1] * n[2] + id[j][1] * n[2] + id[k][2]]
            * (i * xi[0] + (1 - i) * (1 - xi[0]))
            * (j * xi[1] + (1 - j) * (1 - xi[1]))
            * (k * xi[2] + (1 - k) * (1 - xi[2]));
        }
    return v;
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
class UTMField:public Field
{
  int field_id, zone;
  double a, b, n, n2, n3, n4, n5, e, e2, e1, e12, e13, e14, J1, J2, J3, J4,
    Ap, Bp, Cp, Dp, Ep, e4, e6, ep, ep2, ep4, k0, mu_fact;
public:
  UTMField()
  {
    field_id = 1;
    zone = 0;
    options["IField"] = new FieldOptionInt(field_id);
    options["Zone"] = new FieldOptionInt(zone);
    a = 6378137;                /* Equatorial Radius */
    b = 6356752.3142;           /* Rayon Polar Radius */
    /* see http://www.uwgb.edu/dutchs/UsefulData/UTMFormulas.HTM */
    n = (a - b) / (a + b);
    n2 = n * n;
    n3 = n * n * n;
    n4 = n * n * n * n;
    n5 = n * n * n * n * n;
    e = sqrt(1 - b * b / a / a);
    e2 = e * e;
    e1 = (1 - sqrt(1 - e2)) / (1 + sqrt(1 - e2));
    e12 = e1 * e1;
    e13 = e1 * e1 * e1;
    e14 = e1 * e1 * e1 * e1;
    J1 = (3 * e1 / 2 - 27 * e13 / 32);
    J2 = (21 * e12 / 16 - 55 * e14 / 32);
    J3 = 151 * e13 / 96;
    J4 = 1097 * e14 / 512;
    Ap = a * (1 - n + (5. / 4.) * (n2 - n3) + (81. / 64.) * (n4 - n5));
    Bp = -3 * a * n / 2 * (1 - n + (7. / 8.) * (n2 - n3) +
                           (55. / 64.) * (n4 - n5));
    Cp = 14 * a * n2 / 16 * (1 - n + (3. / 4) * (n2 - n3));
    Dp = -35 * a * n3 / 48 * (1 - n + 11. / 16. * (n2 - n3));
    Ep = +315 * a * n4 / 51 * (1 - n);
    e4 = e2 * e2;
    e6 = e2 * e2 * e2;
    ep = e * a / b;
    ep2 = ep * ep;
    ep4 = ep2 * ep2;
    k0 = 0.9996;
    mu_fact = 1 / (k0 * a * (1 - e2 / 4 - 3 * e4 / 64 - 5 * e6 / 256));
  }
  const char *get_name()
  {
    return "UTM";
  }
  double operator() (double x, double y, double z)
  {
    double r = sqrt(x * x + y * y + z * z);
    double lon = atan2(y, x);
    double lat = asin(z / r);
    double meridionalarc = Ap * lat + Bp * sin(2 * lat)
      + Cp * sin(4 * lat) + Dp * sin(6 * lat) + Ep;
    double slat = sin(lat);
    double clat = cos(lat);
    double slat2 = slat * slat;
    double clat2 = clat * clat;
    double clat3 = clat2 * clat;
    double clat4 = clat3 * clat;
    double tlat2 = slat2 / clat2;
    double nu = a / sqrt(1 - e * e * slat2);
    double p = lon - ((zone - 0.5) / 30 - 1) * M_PI;
    double p2 = p * p;
    double p3 = p * p2;
    double p4 = p2 * p2;
    double utm_x =
      k0 * nu * clat * p + (k0 * nu * clat3 / 6) * (1 - tlat2 +
                                                    ep2 * clat2) * p3 + 5e5;
    double utm_y =
      meridionalarc * k0 + k0 * nu * slat * clat / 2 * p2 +
      k0 * nu * slat * clat3 / 24 * (5 - tlat2 + 9 * ep2 * clat2 +
                                     4 * ep4 * clat4) * p4;
    return (*GModel::current()->getFields()->get(field_id)) (utm_x, utm_y, 0);
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
};

class LonLatField:public Field
{
Jean-François Remacle's avatar
Jean-François Remacle committed
  int field_id;
Jean-François Remacle's avatar
Jean-François Remacle committed
  LonLatField()
  {
    field_id = 1;
    options["IField"] = new FieldOptionInt(field_id);
  }
  const char *get_name()
  {
    return "LonLat";
  }
  double operator() (double x, double y, double z)
  {
    return (*GModel::current()->getFields()->get(field_id))
      (atan2(y, x), asin(z / sqrt(x * x + y * y + z * z)), 0);
Jean-François Remacle's avatar
Jean-François Remacle committed
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  double v_in, v_out, x_min, x_max, y_min, y_max, z_min, z_max;
Jean-François Remacle's avatar
Jean-François Remacle committed
  BoxField()
  {
    v_in = v_out = x_min = x_max = y_min = y_max = z_min = z_max = 0;
    options["VIn"] = new FieldOptionDouble(v_in);
    options["VOut"] = new FieldOptionDouble(v_out);
    options["XMin"] = new FieldOptionDouble(x_min);
    options["XMax"] = new FieldOptionDouble(x_max);
    options["YMin"] = new FieldOptionDouble(y_min);
    options["YMax"] = new FieldOptionDouble(y_max);
    options["ZMin"] = new FieldOptionDouble(z_min);
    options["ZMax"] = new FieldOptionDouble(z_max);
  }
  const char *get_name()
  {
    return "Box";
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
  double operator() (double x, double y, double z)
  {
    return (x <= x_max && x >= x_min && y <= y_max && y >= y_min && z <= z_max
            && z >= z_min) ? v_in : v_out;
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  int iField;
  double dmin, dmax, lcmin, lcmax;
Jean-François Remacle's avatar
Jean-François Remacle committed
  const char *get_name()
  {
    return "Threshold";
  }
  ThresholdField()
  {
    iField = 0;
    dmin = 1;
    dmax = 10;
    lcmin = 0.1;
    lcmax = 1;
    options["IField"] = new FieldOptionInt(iField);
    options["DistMin"] = new FieldOptionDouble(dmin);
    options["DistMax"] = new FieldOptionDouble(dmax);
    options["LcMin"] = new FieldOptionDouble(lcmin);
    options["LcMax"] = new FieldOptionDouble(lcmax);
  }
  double operator() (double x, double y, double z)
  {
    Field *field = GModel::current()->getFields()->get(iField);
    double r = ((*field) (x, y, z) - dmin) / (dmax - dmin);
    r = std::max(std::min(r, 1.), 0.);
    double lc = lcmin * (1 - r) + lcmax * r;
    return lc;
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  int iField, kind;
  double delta;
Jean-François Remacle's avatar
Jean-François Remacle committed
  const char *get_name()
  {
    return "Gradient";
  }
  GradientField():iField(0), kind(3), delta(CTX.lc / 1e4)
Jean-François Remacle's avatar
Jean-François Remacle committed
  {
    options["IField"] = new FieldOptionInt(iField);
    options["Kind"] = new FieldOptionInt(kind);
    options["Delta"] = new FieldOptionDouble(delta);
  }
  double operator() (double x, double y, double z)
  {
    Field *field = GModel::current()->getFields()->get(iField);
    double gx, gy, gz;
    switch (kind) {
    case 0:    /* x */
      return ((*field) (x + delta / 2, y, z) -
              (*field) (x - delta / 2, y, z)) / delta;
    case 1:    /* y */
      return ((*field) (x, y + delta / 2, z) -
              (*field) (x, y - delta / 2, z)) / delta;
    case 2:    /* z */
      return ((*field) (x, y, z + delta / 2) -
              (*field) (x, y, z - delta / 2)) / delta;
    case 3:    /* norm */
      gx =
        ((*field) (x + delta / 2, y, z) -
         (*field) (x - delta / 2, y, z)) / delta;
      gy =
        ((*field) (x, y + delta / 2, z) -
         (*field) (x, y - delta / 2, z)) / delta;
      gz =
        ((*field) (x, y, z + delta / 2) -
         (*field) (x, y, z - delta / 2)) / delta;
      return sqrt(gx * gx + gy * gy + gz * gz);
    default:
      Msg::Error("Field %i : Unknown kind (%i) of gradient.", this->id,
Jean-François Remacle's avatar
Jean-François Remacle committed
          kind);
      return MAX_LC;
    }
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
class CurvatureField:public Field
{
  int iField;
  double delta;
public:
  const char *get_name()
  {
    return "Curvature";
  }
  CurvatureField():iField(0),  delta(CTX.lc / 1e4)
  {
    options["IField"] = new FieldOptionInt(iField);
    options["Delta"] = new FieldOptionDouble(delta);
  }
  void grad_norm(Field &f,double x,double y,double z, double *g){
    g[0]=f(x+delta/2,y,z)-f(x-delta/2,y,z);
    g[1]=f(x,y+delta/2,z)-f(x,y-delta/2,z);
    g[2]=f(x,y,z+delta/2)-f(x,y,z-delta/2);
    double n=sqrt(g[0]*g[0]+g[1]*g[1]+g[2]*g[2]);
    g[0]/=n;
    g[1]/=n;
    g[2]/=n;
  }
  double operator() (double x, double y, double z)
  {
    Field *field = GModel::current()->getFields()->get(iField);
    double grad[6][3];
    grad_norm(*field,x+delta/2,y,z,grad[0]);
    grad_norm(*field,x-delta/2,y,z,grad[1]);
    grad_norm(*field,x,y+delta/2,z,grad[2]);
    grad_norm(*field,x,y-delta/2,z,grad[3]);
    grad_norm(*field,x,y,z+delta/2,grad[4]);
    grad_norm(*field,x,y,z-delta/2,grad[5]);
    return (grad[0][0]-grad[1][0]+grad[2][1]-grad[3][1]+grad[4][2]-grad[5][2])/delta;
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
};
Jean-François Remacle's avatar
Jean-François Remacle committed
#if defined(HAVE_GSL)
#include <gsl/gsl_math.h>
#include <gsl/gsl_eigen.h>
class MaxEigenHessianField:public Field
{
  int iField;
  double delta;
  gsl_eigen_symm_workspace *gslwork;
  gsl_matrix *gslmat;
  gsl_vector *eigenvalues;
public:
  const char *get_name()
  {
    return "MaxEigenHessian";
  }
  MaxEigenHessianField():iField(0),  delta(CTX.lc / 1e4)
  {
    options["IField"] = new FieldOptionInt(iField);
    options["Delta"] = new FieldOptionDouble(delta);
    gslwork=gsl_eigen_symm_alloc(3);
    eigenvalues=gsl_vector_alloc(3);
    gslmat=gsl_matrix_alloc(3,3);
  }
  ~MaxEigenHessianField(){
    gsl_eigen_symm_free(gslwork);
    gsl_vector_free(eigenvalues);
    gsl_matrix_free(gslmat);
  }
  double operator() (double x, double y, double z)
  {
    Field *field = GModel::current()->getFields()->get(iField);
    gsl_matrix_set(gslmat,1,0,
          (*field) (x + delta/2 , y+delta/2, z)
        + (*field) (x - delta/2 , y-delta/2, z)
        - (*field) (x - delta/2 , y+delta/2, z)
        - (*field) (x + delta/2 , y-delta/2, z));
    gsl_matrix_set(gslmat,2,0,
          (*field) (x + delta/2 , y, z+delta/2)
        + (*field) (x - delta/2 , y, z-delta/2)
        - (*field) (x - delta/2 , y, z+delta/2)
        - (*field) (x + delta/2 , y, z-delta/2));
    gsl_matrix_set(gslmat,2,1,
          (*field) (x, y + delta/2 , z+delta/2)
        + (*field) (x, y - delta/2 , z-delta/2)
        - (*field) (x, y - delta/2 , z+delta/2)
        - (*field) (x, y + delta/2 , z-delta/2));
    double f=(*field)(x,y,z);
    gsl_matrix_set(gslmat,0,0,
        (*field) (x + delta , y, z)+ (*field) (x - delta , y, z) -2*f);
    gsl_matrix_set(gslmat,1,1,
        (*field) (x, y + delta, z)+ (*field) (x , y - delta, z) -2* f);
    gsl_matrix_set(gslmat,2,2,
        (*field) (x, y ,z + delta)+ (*field) (x , y, z - delta)-2*f);
    gsl_eigen_symm(gslmat,eigenvalues,gslwork);
    return std::max(
      fabs(gsl_vector_get(eigenvalues,0)),
      std::max(
      fabs(gsl_vector_get(eigenvalues,0)),
      fabs(gsl_vector_get(eigenvalues,1))))/(delta*delta);
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
};
#endif

class LaplacianField:public Field
{
  int iField;
  double delta;
public:
  const char *get_name()
  {
    return "Laplacian";
  }
  LaplacianField():iField(0),  delta(CTX.lc / 1e4)
  {
    options["IField"] = new FieldOptionInt(iField);
    options["Delta"] = new FieldOptionDouble(delta);
  }
  double operator() (double x, double y, double z)
  {
    Field *field = GModel::current()->getFields()->get(iField);
      return (
        (*field) (x + delta , y, z)+ (*field) (x - delta , y, z)
        +(*field) (x, y + delta , z)+ (*field) (x, y - delta , z)
        +(*field) (x, y, z + delta )+ (*field) (x, y, z - delta )
        -6* (*field) (x , y, z)
        ) / (delta*delta);
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
};

class MeanField:public Field
{
  int iField;
  double delta;
  int n;
public:
  const char *get_name()
  {
    return "Mean";
  }
  MeanField():iField(0),  delta(CTX.lc / 1e4)
  {
    options["IField"] = new FieldOptionInt(iField);
    options["Delta"] = new FieldOptionDouble(delta);
    //options["N"] = new FieldOptionInt(n);
  }
  double operator() (double x, double y, double z)
  {
    Field *field = GModel::current()->getFields()->get(iField);
    return (
      (*field) (x + delta , y, z)+ (*field) (x - delta , y, z)
      +(*field) (x, y + delta , z)+ (*field) (x, y - delta , z)
      +(*field) (x, y, z + delta )+ (*field) (x, y, z - delta )
      + (*field) (x , y, z)
      ) / 5;
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
};

#if defined(HAVE_MATH_EVAL)
Jean-François Remacle's avatar
Jean-François Remacle committed
  bool error_status;
  std::list < Field * >*list;
  int nvalues;
  char **names;
  double *values;
  void *eval;
Jean-François Remacle's avatar
Jean-François Remacle committed
  int *evaluators_id;
  std::string function;
  char *c_str_function;
Jean-François Remacle's avatar
Jean-François Remacle committed
  double evaluate(double x, double y, double z)
  {
    if(error_status)
      return MAX_LC;
    for(int i = 0; i < nvalues; i++)
    {
      Field *f;
      switch (evaluators_id[i]) {
      case -1:
        values[i] = x;
        break;
        case -2:values[i] = y;
        break;
        case -3:values[i] = z;
        break;
        default:
        {
          f = GModel::current()->getFields()->get(evaluators_id[i]);
          values[i] = f ? (*f) (x, y, z) : MAX_LC;
        }
      }
    }
    return evaluator_evaluate(eval, nvalues, names, values);
  }
  MathEvalExpression() {
    eval = NULL;
    values = NULL;
    c_str_function = NULL;
    evaluators_id = NULL;
  }
  bool set_function(const std::string & f)
  {
    free_members();
    error_status = false;
    c_str_function = strdup(f.c_str());
    eval = evaluator_create(c_str_function);
    if(!eval) {
      error_status = true;
      return false;
    }
    evaluator_get_variables(eval, &names, &nvalues);
    values = new double[nvalues];
    evaluators_id = new int[nvalues];
    for(int i = 0; i < nvalues; i++) {
      int id;
      if(!strcmp("x", names[i]))
        evaluators_id[i] = -1;
      else if(!strcmp("y", names[i]))
        evaluators_id[i] = -2;
      else if(!strcmp("z", names[i]))
        evaluators_id[i] = -3;
      else if(sscanf(names[i], "F%i", &id) == 1)
        evaluators_id[i] = id;
      else {
        Msg::Error("Unknown matheval argument \"%s\"\n", names[i]);
Jean-François Remacle's avatar
Jean-François Remacle committed
        error_status = true;
        return false;
      }
    }
    return true;
  }
  void free_members()
  {
    if(c_str_function)
      free(c_str_function);
    if(eval)
      evaluator_destroy(eval);
    if(values)
      delete[]values;
    if(evaluators_id)
      delete evaluators_id;
  }
  ~MathEvalExpression() {
    free_members();
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  MathEvalExpression expr;
  std::string f;
Jean-François Remacle's avatar
Jean-François Remacle committed
  MathEvalField()
  {
    options["F"] = new FieldOptionString(f, &update_needed);
  }
  double operator() (double x, double y, double z)
  {
    if(update_needed) {
      if(!expr.set_function(f))
        Msg::Error("Field %i : Invalid matheval expression \"%s\"\n",
Jean-François Remacle's avatar
Jean-François Remacle committed
            this->id, f.c_str());
      update_needed = false;
    }
    return expr.evaluate(x, y, z);
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
  const char *get_name()
  {
    return "MathEval";
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  MathEvalExpression expr[3];
  std::string f[3];
  int ifield;
Jean-François Remacle's avatar
Jean-François Remacle committed
  ParametricField()
  {
    options["IField"] = new FieldOptionInt(ifield);
    options["FX"] = new FieldOptionString(f[0], &update_needed);
    options["FY"] = new FieldOptionString(f[1], &update_needed);
    options["FZ"] = new FieldOptionString(f[2], &update_needed);
  }
  double operator() (double x, double y, double z)
  {
    if(update_needed) {
      for(int i = 0; i < 3; i++) {
        if(!expr[i].set_function(f[i]))
          Msg::Error("Field %i : Invalid matheval expression \"%s\"\n",
Jean-François Remacle's avatar
Jean-François Remacle committed
              this->id, f[i].c_str());
      }
      update_needed = false;
    }
    return (*GModel::current()->getFields()->get(ifield))
Christophe Geuzaine's avatar
Christophe Geuzaine committed
      (expr[0].evaluate(x, y, z), expr[1].evaluate(x, y, z),
       expr[2].evaluate(x, y, z));
Jean-François Remacle's avatar
Jean-François Remacle committed
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
  const char *get_name()
  {
    return "Param";
  }
  OctreePost *octree;
Jean-François Remacle's avatar
Jean-François Remacle committed
  double operator() (double x, double y, double z)
  {
    // FIXME: should test unique view num instead, but that would be slower
Jean-François Remacle's avatar
Jean-François Remacle committed
    if(view_index < 0 || view_index >= (int)PView::list.size())
      return MAX_LC;
Jean-François Remacle's avatar
Jean-François Remacle committed
      if(octree)
        delete octree;
      octree = new OctreePost(PView::list[view_index]);
      update_needed = false;
    }
    double l = 0.;
    if(!octree->searchScalar(x, y, z, &l, 0)) {
      // uncomment the following to try really hard to find an element
      // around the point
      /*
         double fact[9] = {0.001, 0.0025, 0.005, 0.0075, 0.01, 0.025, 0.05, 0.075, 0.1};
         for(int i = 0; i < 9; i++){
         double eps = CTX.lc * fact[i];
         if(octree->searchScalar(x + eps, y, z, &l, 0)) break;
         if(octree->searchScalar(x - eps, y, z, &l, 0)) break;
         if(octree->searchScalar(x, y + eps, z, &l, 0)) break;
         if(octree->searchScalar(x, y - eps, z, &l, 0)) break;
         if(octree->searchScalar(x, y, z + eps, &l, 0)) break;
         if(octree->searchScalar(x, y, z - eps, &l, 0)) break;
         if(octree->searchScalar(x + eps, y - eps, z - eps, &l, 0)) break;
         if(octree->searchScalar(x + eps, y + eps, z - eps, &l, 0)) break;
         if(octree->searchScalar(x - eps, y - eps, z - eps, &l, 0)) break;
         if(octree->searchScalar(x - eps, y + eps, z - eps, &l, 0)) break;
         if(octree->searchScalar(x + eps, y - eps, z + eps, &l, 0)) break;
         if(octree->searchScalar(x + eps, y + eps, z + eps, &l, 0)) break;
         if(octree->searchScalar(x - eps, y - eps, z + eps, &l, 0)) break;
         if(octree->searchScalar(x - eps, y + eps, z + eps, &l, 0)) break;
         }
       */
    }
    if(l <= 0) return MAX_LC;
Jean-François Remacle's avatar
Jean-François Remacle committed
    return l;
  }
  const char *get_name()
  {
    return "PostView";
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
    octree = NULL;
    options["IView"] = new FieldOptionInt(view_index, &update_needed);
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
    if(octree)
      delete octree;
  }
  FieldDialogBox *&dialog_box()
  {
    static FieldDialogBox *dialogBox = NULL;
    return dialogBox;
  }
Jean-François Remacle's avatar
Jean-François Remacle committed
  std::list < int >idlist;
Jean-François Remacle's avatar
Jean-François Remacle committed
  MinField()
  {
    options["FieldsList"] = new FieldOptionList(idlist, &update_needed);
  }
  double operator() (double x, double y, double z)
  {
    double v = MAX_LC;
    for(std::list < int >::iterator it = idlist.begin(); it != idlist.end();
        it++) {
      Field *f = (GModel::current()->getFields()->get(*it));
      if(f)
        v = std::min(v, (*f) (x, y, z));
    }