Hi Romain,

please have a look at the compiler messages:

tutorialproblem_coupled.hh:173:14: erreur: déclaration ne déclarant rien du tout tutorialproblem_coupled.hh:173:33: erreur: 'fvElemGeom' was not declared in this scope tutorialproblem_coupled.hh:174:14: erreur: 'VolumeVariables' was not declared in this scope tutorialproblem_coupled.hh:174:31: erreur: 'volVars' was not declared in this scope
The first semicolon in each of the lines
FVElementGeometry; fvElemGeom;
             VolumeVariables; volVars;
is too much.

tutorialproblem_coupled.hh:176:14: erreur: 'ElementIterator' was not declared in this scope
tutorialproblem_coupled.hh:176:30: erreur: expected ';' before 'elemIt'
tutorialproblem_coupled.hh:177:30: erreur: expected ';' before 'elemEndIt'
tutorialproblem_coupled.hh:178:21: erreur: 'elemIt' was not declared in this scope tutorialproblem_coupled.hh:178:31: erreur: 'elemEndIt' was not declared in this scope
You have to declare the type "ElementIterator." You can copy the typedef from heterogeneousproblem.hh.

tutorialproblem_coupled.hh:186:38: erreur: 'numVerts' was not declared in this scope
You have to uncomment the line defining numVerts (184?).

tutorialproblem_coupled.hh:198:22: erreur: 'Kzz' was not declared in this scope
You have to comment this line since you only want to provide Kxx.

Kind regards
Bernd

..
1- But I don't how to declare variables and where to do it. Because in the " "test/boxmodels/co2/heterogeneousproblem.hh.", variables are only declared in the routine and it compiles but not when I copy it in "tutorialproblem_coupled.hh". 2- If I understood well it is from intrinsicPermeability that my permeability will come from "tutorialspatialparams.hh" ?
Thanks to the community,
Cheers,
Romain.
ps : the addOutputVtkFields :
void addOutputVtkFields()
        {
typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;
            unsigned numVertices = this->gridView().size(dim);


//             //create required scalar fields
ScalarField *Kxx = this->resultWriter().allocateManagedBuffer(numVertices); // ScalarField *boxVolume = this->resultWriter().allocateManagedBuffer(numVertices);
// //
// //            (*boxVolume) = 0;
// //
            //Fill the scalar fields with values
              unsigned numElements = this->gridView().size(0);
ScalarField *rank = this->resultWriter().allocateManagedBuffer(numElements);
//
             FVElementGeometry; fvElemGeom;
             VolumeVariables; volVars;
ElementIterator elemIt = this->gridView().template begin<0>();
ElementIterator elemEndIt = this->gridView().template end<0>();
             for (; elemIt != elemEndIt; ++elemIt)
             {
                 int idx = this->elementMapper().map(*elemIt);
                 (*rank)[idx] = this->gridView().comm().rank();
                 fvElemGeom.update(this->gridView(), *elemIt);
//
//                 int numVerts = elemIt->template count<dim> ();
//
                 for (int i = 0; i < numVerts; ++i)
                 {
int globalIdx = this->vertexMapper().map(*elemIt, i, dim);
volVars.update(this->model().curSol()[globalIdx],
                                    *this,
                                    *elemIt,
                                    fvElemGeom,
                                    i,
                                    false);
// (*boxVolume)[globalIdx] += fvElemGeom.subContVol[i].volume; Scalar perm = this->spatialParams().intrinsicPermeability(*elemIt, fvElemGeom, i);
                    (*Kxx)[globalIdx] = perm;
                   (*Kzz)[globalIdx] =   perm[dim-1][dim-1];
                 }
             }
//pass the scalar fields to the vtkwriter
// this->resultWriter().attachVertexData(*boxVolume, "boxVolume");
            this->resultWriter().attachVertexData(*Kxx, "Kxx");
    };

__________________________

/Avant d'imprimer, pensez à l'environnement ! Please consider the environment before printing ! / /Ce message et toutes ses pièces jointes sont confidentiels et établis à l'intention exclusive de ses destinataires. Toute utilisation non conforme à sa destination, toute diffusion ou toute publication, totale ou partielle, est interdite, sauf autorisation expresse. IFP Energies nouvelles décline toute responsabilité au titre de ce message. This message and any attachments are confidential and intended solely for the addressees. Any unauthorised use or dissemination is prohibited. IFP Energies nouvelles should not be liable for this message./
__________________________



_______________________________________________
Dumux mailing list
Dumux@listserv.uni-stuttgart.de
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux


--
_____________________________________________________________________

Bernd Flemisch                               phone: +49 711 685 69162
IWS, Universität Stuttgart                   fax:   +49 711 685 60430
Pfaffenwaldring 61                  email: be...@iws.uni-stuttgart.de
D-70569 Stuttgart                  url: www.hydrosys.uni-stuttgart.de
_____________________________________________________________________

_______________________________________________
Dumux mailing list
Dumux@listserv.uni-stuttgart.de
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux

Reply via email to