Dear DuMuX developer,

I try to write gas-phase saturations of all nodes (in 2p2c model) in a *.txt 
file for each time integration, and  use them next in another code to conduct 
some geochemical simulations.
I followed the way that results are written in vtk files and prepared  void 
postTimeStep as shown below (I didn't show it here, but (*Saturation_...) are 
all written into a text file).
It works perfectly when I run my simulations on one (1) processor. However, on 
several processors, I obtain a larger global Idx than those I find in vtk files 
for each rank. For example, in rank 0, maximum global Idx is 7340 in *.txt file 
while in vtk file is only  5251. I tried to figure out what is the cause of 
this difference, and I mentioned that in *.txt files for all nodes located on 
the interface of two ranks (processors), I have two different global Idxs, for 
example (I give here the value of saturation on a vertex:
globalIdx, x(m), y(m), z(m), and gas-phase saturation(-):
445,  -10.53,  24,  4.8,  0.55
and
5346,  -10.53,  24,  4.8,  0

As can be seen, for the same global position I have two different Idxs. How can 
it be possible? How can we have two different indxs for the same global 
position? The gas-phase saturation of  globalIdx 445 corresponds to the value 
stored in the vtk file for the same position (-10.53,  24,  4.8). What does the 
second index represent? Is there any problem in the code that I used for 
writing the results?

Hope that I was clear in my explanations
Please let me know if you need more details
Thanks in advance for your help
Best regard
Ali NOWAMOOZ

void postTimeStep() {

FVElementGeometry fvGeometry;
ElementVolumeVariables elemVolVars;
// get the number of degrees of freedom
unsigned numDofs = this->model().numDofs();
typedef Dune::BlockVector<Dune::FieldVector<double, 1> > ScalarField;
ScalarField *Saturation_nw = 
this->resultWriter().allocateManagedBuffer(numDofs);
ScalarField *Saturation_x = this->resultWriter().allocateManagedBuffer(numDofs);
ScalarField *Saturation_y = this->resultWriter().allocateManagedBuffer(numDofs);
ScalarField *Saturation_z = this->resultWriter().allocateManagedBuffer(numDofs);
ScalarField *Saturation_scv = 
this->resultWriter().allocateManagedBuffer(numDofs);
(*Saturation_nw) = 0;
(*Saturation_x) = 0;
(*Saturation_y) = 0;
(*Saturation_z) = 0;
(*Saturation_scv) = 0;

ElementIterator eIt = this->gridView().template begin<0>();
ElementIterator eEndIt = this->gridView().template end<0>();

for (; eIt != eEndIt; ++eIt){
elemVolVars.update(*this,*eIt,fvGeometry,false);
fvGeometry.update(this->gridView(), *eIt);

for (int scvIdx = 0; scvIdx < fvGeometry.numScv; ++scvIdx){
int gIobaldx = this->model().dofMapper().map(*eIt, scvIdx, dofCodim);
const GlobalPosition gpos = eIt->geometry().corner(scvIdx);
(*Saturation_nw)[ gIobaldx] = elemVolVars[scvIdx].saturation(gPhaseIdx);
(*Saturation_x)[ gIobaldx] = gpos[0];
(*Saturation_y)[ gIobaldx] = gpos[1];
(*Saturation_z)[ gIobaldx] = gpos[2];
(*Saturation_scv)[ gIobaldx] += 1;
}

}

}

_______________________________________________
Dumux mailing list
[email protected]
https://listserv.uni-stuttgart.de/mailman/listinfo/dumux

Reply via email to