SolveSuperLU (const MatriceMorse<R> &AA, int strategy, double ttgv, double epsilon,
		              double pivot, double pivot_sym, string &param_char, KN<long> pperm_r,
		              KN<long> pperm_c):
			eps(epsilon), epsr(0),
			tgv(ttgv),
			etree(0), string_option(param_char), perm_r(pperm_r), perm_c(pperm_c),
			RR(0), CC(0),
			tol_pivot_sym(pivot_sym), tol_pivot(pivot) {
			SuperMatrix B, X;
			SuperLUStat_t stat;
			void *work = 0;
			int info, lwork = 0/*, nrhs = 1*/;
			int i;
			double ferr[1];
			double berr[1];
			double rpg, rcond;
			R *bb;
			R *xx;

			A.Store = 0;
			B.Store = 0;
			X.Store = 0;
			L.Store = 0;
			U.Store = 0;

			int status;

			n = AA.n;
			m = AA.m;
			nnz = AA.nbcoef;

			arow = AA.a;
			asubrow = AA.cl;
			xarow = AA.lg;

			/* FreeFem++ use Morse Format */
			// FFCS - "this->" required by g++ 4.7
			this->CompRow_to_CompCol(m, n, nnz, arow, asubrow, xarow,
			                         &a, &asub, &xa);

			/* Defaults */
			lwork = 0;
			// nrhs = 0;

			/* Set the default values for options argument:
			 *  options.Fact = DOFACT;
			 *  options.Equil = YES;
			 *  options.ColPerm = COLAMD;
			 *  options.DiagPivotThresh = 1.0;
			 *  options.Trans = NOTRANS;
			 *  options.IterRefine = NOREFINE;
			 *  options.SymmetricMode = NO;
			 *  options.PivotGrowth = NO;
			 *  options.ConditionNumber = NO;
			 *  options.PrintStat = YES;
			 */
			set_default_options(&options);

			printf(".. default options:\n");
			printf("\tFact\t %8d\n", options.Fact);
			printf("\tEquil\t %8d\n", options.Equil);
			printf("\tColPerm\t %8d\n", options.ColPerm);
			printf("\tDiagPivotThresh %8.4f\n", options.DiagPivotThresh);
			printf("\tTrans\t %8d\n", options.Trans);
			printf("\tIterRefine\t%4d\n", options.IterRefine);
			printf("\tSymmetricMode\t%4d\n", options.SymmetricMode);
			printf("\tPivotGrowth\t%4d\n", options.PivotGrowth);
			printf("\tConditionNumber\t%4d\n", options.ConditionNumber);
			printf("..\n");

			if (!string_option.empty()) {read_options_freefem(string_option, &options);}

			printf(".. options:\n");
			printf("\tFact\t %8d\n", options.Fact);
			printf("\tEquil\t %8d\n", options.Equil);
			printf("\tColPerm\t %8d\n", options.ColPerm);
			printf("\tDiagPivotThresh %8.4f\n", options.DiagPivotThresh);
			printf("\tTrans\t %8d\n", options.Trans);
			printf("\tIterRefine\t%4d\n", options.IterRefine);
			printf("\tSymmetricMode\t%4d\n", options.SymmetricMode);
			printf("\tPivotGrowth\t%4d\n", options.PivotGrowth);
			printf("\tConditionNumber\t%4d\n", options.ConditionNumber);
			printf("..\n");

			Dtype_t R_SLU = SuperLUDriver<R>::R_SLU_T();

			// FFCS - "this->" required by g++ 4.7
			this->Create_CompCol_Matrix(&A, m, n, nnz, a, asub, xa, SLU_NC, R_SLU, SLU_GE);

			this->Create_Dense_Matrix(&B, m, 0, (R *)0, m, SLU_DN, R_SLU, SLU_GE);
			this->Create_Dense_Matrix(&X, m, 0, (R *)0, m, SLU_DN, R_SLU, SLU_GE);

			if (etree.size() == 0) {etree.resize(n);}

			if (perm_r.size() == 0) {perm_r.resize(n);}

			if (perm_c.size() == 0) {perm_c.resize(n);}

			if (!(RR = new double[n])) {
				ABORT("SUPERLU_MALLOC fails for R[].");
			}

			for (int ii = 0; ii < n; ii++) {
				RR[ii] = 1.;
			}

			if (!(CC = new double[m])) {
				ABORT("SUPERLU_MALLOC fails for C[].");
			}

			for (int ii = 0; ii < n; ii++) {
				CC[ii] = 1.;
			}

			ferr[0] = 0;
			berr[0] = 0;
			/* Initialize the statistics variables. */
			StatInit(&stat);

			/* ONLY PERFORM THE LU DECOMPOSITION */
			B.ncol = 0;	/* Indicate not to solve the system */
			SuperLUDriver<R>::gssvx(&options, &A, perm_c, perm_r, etree, equed, RR, CC,
			                        &L, &U, work, lwork, &B, &X, &rpg, &rcond, ferr, berr, &Glu,
			                        &mem_usage, &stat, &info);

			if (verbosity > 2) {
				printf("LU factorization: dgssvx() returns info %d\n", info);
			}

			if (verbosity > 3) {
				if (info == 0 || info == n + 1) {
					if (options.PivotGrowth) {printf("Recip. pivot growth = %e\n", rpg);}

					if (options.ConditionNumber) {
						printf("Recip. condition number = %e\n", rcond);
					}

					Lstore = (SCformat *)L.Store;
					Ustore = (NCformat *)U.Store;
					printf("No of nonzeros in factor L = %d\n", Lstore->nnz);
					printf("No of nonzeros in factor U = %d\n", Ustore->nnz);
					printf("No of nonzeros in L+U = %d\n", Lstore->nnz + Ustore->nnz - n);
					printf("L\\U MB %.3f\ttotal MB needed %.3f\texpansions %d\n",
					       mem_usage.for_lu / 1e6, mem_usage.total_needed / 1e6,
					       stat.expansions
					       );
					fflush(stdout);
				} else if (info > 0 && lwork == -1) {
					printf("** Estimated memory: %d bytes\n", info - n);
				}
			}

			if (verbosity > 5) {StatPrint(&stat);}

			StatFree(&stat);
			if (B.Store) {Destroy_SuperMatrix_Store(&B);}

			if (X.Store) {Destroy_SuperMatrix_Store(&X);}

			options.Fact = FACTORED;/* Indicate the factored form of A is supplied. */
		}
AnyType yams_Op::operator () (Stack stack)  const {
	// initialisation
	MeshPoint *mp(MeshPointStack(stack)), mps = *mp;
	Mesh3 *pTh = GetAny<Mesh3 *>((*eTh)(stack));

	ffassert(pTh);
	Mesh3 &Th3 = *pTh;
	int nv = Th3.nv;
	int nt = Th3.nt;
	int nbe = Th3.nbe;

	KN<int> defaultintopt(23);
	KN<double> defaultfopt(14);
	defaultintopt = 0;
	defaultfopt = 0.;
	yams_inival(defaultintopt, defaultfopt);

	KN<int> intopt(23);

	for (int ii = 0; ii < 23; ii++) {
		intopt[ii] = defaultintopt[ii];
	}

	KN<double> fopt(14);

	for (int ii = 0; ii < 14; ii++) {
		fopt[ii] = defaultfopt[ii];
	}

	assert(fopt.N() == 14);

	if (nargs[0]) {
		KN<int> intopttmp = GetAny<KN_<long> >((*nargs[0])(stack));
		if (intopttmp.N() != 13) {
			cerr << "the size of vector loptions is 13 " << endl;
			exit(1);
		} else {
			for (int ii = 0; ii < 13; ii++) {
				intopt[wrapper_intopt[ii]] = intopttmp[ii];
			}
		}
	}

	if (nargs[1]) {
		KN<double> fopttmp = GetAny<KN_<double> >((*nargs[1])(stack));
		if (fopttmp.N() != 11) {
			cerr << "the size of vector loptions is 11 not " << fopttmp.N() << endl;
			ExecError("FreeYams");
		} else {
			for (int ii = 0; ii < 11; ii++) {
				fopt[wrapper_fopt[ii]] = fopttmp[ii];
			}
		}
	}

	intopt[0] = arg(3, stack, intopt[0] != 1);
	intopt[8] = arg(4, stack, intopt[8]);
	fopt[7] = arg(5, stack, fopt[7]);
	fopt[8] = arg(6, stack, fopt[7]);
	fopt[6] = arg(7, stack, fopt[6]);
	intopt[22] = arg(8, stack, intopt[22]);	// optim option
	if (nargs[9]) {intopt[17] = 1;}

	fopt[13] = arg(9, stack, fopt[13]);	// ridge angle
	intopt[21] = arg(10, stack, intopt[21]);// absolue
	intopt[11] = arg(11, stack, (int)verbosity);// verbosity
	intopt[17] = arg(12, stack, intopt[17]);// no ridge
	intopt[18] = arg(13, stack, intopt[18]);// nb smooth
	if (verbosity > 1) {
		cout << " fopt = [";

		for (int i = 0; i < 11; ++i) {
			cout << fopt[wrapper_fopt[i]] << (i < 10 ? "," : "];\n");
		}

		cout << " intopt = [";

		for (int i = 0; i < 13; ++i) {
			cout << intopt[wrapper_intopt[i]] << (i < 12 ? "," : "];\n");
		}
	}

	/*
	 * KN<int> intopt(arg(0,stack,defaultintopt));
	 * assert( intopt.N() == 23 );
	 * KN<double> fopt(arg(1,stack,defaultfopt));
	 * assert( fopt.N() == 14 );
	 */
	KN<double> metric;

	int mtype = type;
	if (nargs[2]) {
		metric = GetAny<KN_<double> >((*nargs[2])(stack));
		if (metric.N() == Th3.nv) {
			mtype = 1;
			intopt[1] = 0;
		} else if (metric.N() == 6 * Th3.nv) {
			intopt[1] = 1;
			mtype = 3;
		} else {
			cerr << "sizeof vector metric is incorrect, size will be Th.nv or 6*Th.nv" << endl;
		}
	} else if (nbsol > 0) {
		if (type == 1) {
			intopt[1] = 0;
			metric.resize(Th3.nv);
			metric = 0.;
		} else if (type == 3) {
			intopt[1] = 1;
			metric.resize(6 * Th3.nv);
			metric = 0.;
		}
	} else {
		if (intopt[1] == 0) {metric.resize(Th3.nv); metric = 0.;} else if (intopt[1] == 1) {metric.resize(6 * Th3.nv); metric = 0.;}
	}

	// mesh for yams
	yams_pSurfMesh yamsmesh;
	yamsmesh = (yams_pSurfMesh)calloc(1, sizeof(yams_SurfMesh));
	if (!yamsmesh) {
		cerr << "allocation error for SurfMesh for yams" << endl;
	}

	yamsmesh->infile = NULL;
	yamsmesh->outfile = NULL;
	yamsmesh->type = M_SMOOTH | M_QUERY | M_DETECT | M_BINARY | M_OUTPUT;

	mesh3_to_yams_pSurfMesh(Th3, intopt[8], intopt[22], yamsmesh);

	// solution for freeyams2
	if (nbsol) {
		MeshPoint *mp3(MeshPointStack(stack));

		KN<bool> takemesh(nv);
		takemesh = false;

		for (int it = 0; it < nt; it++) {
			for (int iv = 0; iv < 4; iv++) {
				int i = Th3(it, iv);

				if (takemesh[i] == false) {
					mp3->setP(&Th3, it, iv);

					for (int ii = 0; ii < nbsolsize; ii++) {
						metric[i * nbsolsize + ii] = GetAny<double>((*sol[ii])(stack));
					}

					takemesh[i] = true;
				}
			}
		}
	}

	if (verbosity > 10) {
		cout << "nbsol  " << nargs[2] << endl;
	}

	if (nargs[2] || (nbsol > 0)) {
		float hmin, hmax;
		solyams_pSurfMesh(yamsmesh, mtype, metric, hmin, hmax);
		yamsmesh->nmfixe = yamsmesh->npfixe;
		if (fopt[7] < 0.0) {
			fopt[7] = max(fopt[7], hmin);
		}

		if (fopt[8] < 0.0) {
			fopt[8] = max(fopt[8], hmax);
		}
	} else {
		yamsmesh->nmfixe = 0;
	}

	int infondang = 0, infocc = 0;
	int res = yams_main(yamsmesh, intopt, fopt, infondang, infocc);
	if (verbosity > 10) {
		cout << " yamsmesh->dim " << yamsmesh->dim << endl;
	}

	if (res > 0) {
		cout << " problem with yams :: error " << res << endl;
		ExecError("Freeyams error");
	}

	Mesh3 *Th3_T = yams_pSurfMesh_to_mesh3(yamsmesh, infondang, infocc, intopt[22]);

	// recuperer la solution ????
	if (verbosity > 10) {
		cout << &yamsmesh->point << " " << &yamsmesh->tria << " " << &yamsmesh->geom << " " << &yamsmesh->tgte << endl;
		cout << &yamsmesh << endl;
	}

	free(yamsmesh->point);
	free(yamsmesh->tria);
	free(yamsmesh->geom);
	free(yamsmesh->tgte);
	if (yamsmesh->metric) {free(yamsmesh->metric);}

	if (yamsmesh->edge) {free(yamsmesh->edge);}

	if (yamsmesh->tetra) {free(yamsmesh->tetra);}

	free(yamsmesh);

	*mp = mps;
	Add2StackOfPtr2FreeRC(stack, Th3_T);
	return SetAny<pmesh3>(Th3_T);
}