Skip to content
Open
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
13 changes: 11 additions & 2 deletions femtools/VTK_interfaces.F90
Original file line number Diff line number Diff line change
Expand Up @@ -103,7 +103,7 @@ end subroutine vtkwriteghostlevels

contains

subroutine vtk_write_state(filename, index, model, state, write_region_ids, write_columns, stat)
subroutine vtk_write_state(filename, index, model, state, write_region_ids, write_columns, metadata, stat)
!!< Write the state variables out to a vtu file. Two different elements
!!< are supported along with fields corresponding to each of them.
!!<
Expand All @@ -118,6 +118,7 @@ subroutine vtk_write_state(filename, index, model, state, write_region_ids, writ
type(state_type), dimension(:), intent(in) :: state
logical, intent(in), optional :: write_region_ids
logical, intent(in), optional :: write_columns
type(scalar_field), dimension(:), optional :: metadata
integer, intent(out), optional :: stat
type(mesh_type), pointer :: model_mesh

Expand Down Expand Up @@ -214,12 +215,14 @@ subroutine vtk_write_state(filename, index, model, state, write_region_ids, writ
tfields=ltfields, &
write_region_ids=write_region_ids, &
write_columns=write_columns, &
metadata=metadata, &
stat=stat)

end subroutine vtk_write_state

subroutine vtk_write_fields(filename, index, position, model, sfields,&
& vfields, tfields, write_region_ids, write_columns, number_of_partitions, stat)
& vfields, tfields, write_region_ids, write_columns,&
& metadata, number_of_partitions, stat)
!!< Write the state variables out to a vtu file. Two different elements
!!< are supported along with fields corresponding to each of them.
!!<
Expand All @@ -238,6 +241,7 @@ subroutine vtk_write_fields(filename, index, position, model, sfields,&
type(tensor_field), dimension(:), intent(in), optional :: tfields
logical, intent(in), optional :: write_region_ids
logical, intent(in), optional :: write_columns
type(scalar_field), dimension(:), optional :: metadata
Copy link
Contributor

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Hm, this seems to abuse scalar_field as a container for a real array (with a name). A scalar_field is defined on a mesh (function space!), and it suggests we care which mesh. I can see this is a compact way of supplying one or more real arrays of different lengths with names, but if people want to use this in a generic way they really will have to abuse scalar_fields by allocating a %val array that is inconsistent with whatever random mesh they choose to use. Also when first reading this code, I thought it was just writing out an extra field (i.e an actual spatial field), so I think it's better to make it explcit that we're just supplying a real array of arbitrary length.

Since we only use this for time at the moment, I would be happy to just have a real, intent(in):: time_value argument (where the name is hardcoded to TimeValue). Alternatively you could use an array of real_vectors (see femtools/Futils.F90) with a separate array of names.

!!< If present, only write for processes 1:number_of_partitions (assumes the other partitions are empty)
integer, optional, intent(in):: number_of_partitions
integer, intent(out), optional :: stat
Expand Down Expand Up @@ -758,6 +762,11 @@ subroutine vtk_write_fields(filename, index, position, model, sfields,&

end if

if (present(metadata)) then
do i=1, size(metadata)
call vtkwritefd(metadata(i)%val, trim(metadata(i)%name))
end do
end if

!----------------------------------------------------------------------
! Close the file
Expand Down
17 changes: 15 additions & 2 deletions femtools/Write_State.F90
Original file line number Diff line number Diff line change
Expand Up @@ -27,7 +27,7 @@ module write_state_module
vtk_write_state_new_options

! Static variables set by update_dump_times and used by do_write_state
logical, save :: last_times_initialised = .false.
logical, save :: last_times_initialised = .false., time_as_metadata = .false.
real, save :: last_dump_time
real, save :: last_dump_cpu_time
real, save :: last_dump_wall_time
Expand All @@ -41,6 +41,8 @@ subroutine initialise_write_state
!!< variables)

call update_dump_times

time_as_metadata=have_option('/io/time_as_metadata')

end subroutine initialise_write_state

Expand Down Expand Up @@ -241,6 +243,7 @@ subroutine vtk_write_state_new_options(filename, index, state, write_region_ids)
type(scalar_field), dimension(:), allocatable :: lsfields
type(vector_field), dimension(:), allocatable :: lvfields
type(tensor_field), dimension(:), allocatable :: ltfields
type(scalar_field), dimension(:), allocatable :: metadata
character(len = FIELD_NAME_LEN) :: field_name, mesh_name
integer :: i, f, counter
logical :: multi_state
Expand Down Expand Up @@ -351,14 +354,24 @@ subroutine vtk_write_state_new_options(filename, index, state, write_region_ids)
ewrite(2, "(a,i0,a)") "Writing ", size(ltfields), " tensor field(s)"

model_coordinate=>get_external_coordinate_field(state(1), model_mesh)

if (time_as_metadata) then
allocate(metadata(1))
metadata(1)=extract_scalar_field(state(1),"Time")
metadata%name = "TimeValue"
else
allocate(metadata(0))
end if

call vtk_write_fields(filename, index, &
model_coordinate, &
model_mesh, &
sfields=lsfields, &
vfields=lvfields, &
tfields=ltfields, &
write_region_ids=write_region_ids)
write_region_ids=write_region_ids, &
metadata = metadata &
)

ewrite(1, *) "Exiting vtk_write_state_new_options"

Expand Down
75 changes: 72 additions & 3 deletions libvtkfortran/fvtkfortran.F90
Original file line number Diff line number Diff line change
Expand Up @@ -56,7 +56,7 @@ module vtkfortran

public :: vtkopen, vtkclose, vtkpclose, vtkwritemesh, vtkwritesn,&
& vtkwritesc, vtkwritevn, vtkwritevc, vtkwritetn, vtkwritetc, &
& vtksetactivescalars, vtksetactivevectors, &
& vtkwritefd, vtksetactivescalars, vtksetactivevectors, &
& vtksetactivetensors

interface vtkopen
Expand Down Expand Up @@ -287,7 +287,45 @@ SUBROUTINE VTKWRITEDTC_C(v1, v2, v3, v4, v5, v6, v7, v8, v9, name,&
end SUBROUTINE VTKWRITEDTC_C
module procedure vtkwriteftc_f90, vtkwritedtc_f90
end interface


interface vtkwritefd
! Write generic field (field-based) to the current vtk file.
SUBROUTINE VTKWRITEFFD_C(val, vlen, name,&
& len) bind(c,name="vtkwriteffd")
use iso_c_binding
implicit none
REAL(c_float) :: val(*)
character(kind=c_char,len=1), dimension(*) :: name
integer(kind=c_int) :: vlen, len
end SUBROUTINE VTKWRITEFFD_C
SUBROUTINE VTKWRITEDFD_C(val, vlen, name,&
& len) bind(c,name="vtkwritedfd")
use iso_c_binding
implicit none
REAL(c_double) :: val(*)
character(kind=c_char,len=1), dimension(*) :: name
integer(kind=c_int) :: vlen, len
end SUBROUTINE VTKWRITEDFD_C
SUBROUTINE VTKWRITEIFD_C(val, vlen, name,&
& len) bind(c,name="vtkwriteifd")
use iso_c_binding
implicit none
integer(c_int) :: val(*)
character(kind=c_char,len=1), dimension(*) :: name
integer(kind=c_int) :: vlen, len
end SUBROUTINE VTKWRITEIFD_C
SUBROUTINE VTKWRITECFD_C(val, vlen, name,&
& len) bind(c,name="vtkwritecfd")
use iso_c_binding
implicit none
character(kind=c_char,len=1), dimension(*) :: val, name
integer(kind=c_int) :: vlen, len
end SUBROUTINE VTKWRITECFD_C

module procedure vtkwriteffd_f90, vtkwritedfd_f90, &
vtkwriteifd_f90, vtkwritecfd_f90
end interface

interface vtksetactivescalars
subroutine vtksetactivescalars_c(name, len) bind(c,name="vtksetactivescalars")
use iso_c_binding
Expand Down Expand Up @@ -443,7 +481,38 @@ subroutine vtkwritedtc_f90(v1, v2, v3, v4, v5, v6, v7, v8, v9, name)

call vtkwritedtc_c(v1, v2, v3, v4, v5, v6, v7, v8, v9, name, len(name))

end subroutine vtkwritedtc_f90
end subroutine vtkwritedtc_f90

subroutine vtkwriteifd_f90(val, name)
integer(c_int) :: val(:)
character(len=*) :: name

call vtkwriteifd_c(val, size(val), name, len(name))

end subroutine vtkwriteifd_f90

subroutine vtkwriteffd_f90(val, name)
real(c_float) :: val(:)
character(len=*) :: name

call vtkwriteffd_c(val, size(val), name, len(name))

end subroutine vtkwriteffd_f90

subroutine vtkwritedfd_f90(val, name)
real(c_double) :: val(:)
character(len=*) :: name

call vtkwritedfd_c(val, size(val), name, len(name))

end subroutine vtkwritedfd_f90

subroutine vtkwritecfd_f90(val, name)
character(len=*) :: val, name

call vtkwritecfd_c(val, len(val), name, len(name))

end subroutine vtkwritecfd_f90

subroutine vtksetactivescalars_f90(name)
character(len=*) name
Expand Down
94 changes: 94 additions & 0 deletions libvtkfortran/vtkfortran.cpp
Original file line number Diff line number Diff line change
Expand Up @@ -655,6 +655,100 @@ extern "C" {
return;
}

/**
Writes integer field data values.
@param[in] vect Data array.
@param[in] vlen length of vect.
@param[in] name Variable name to be written to metadata.
@param[in] len Length of variable name.
*/
void vtkwriteifd(int *vect, int *vlen, char *name, int *len){
string tag(name, *len);
vtkIntArray *newField = vtkIntArray::New();
newField->SetName( tag.c_str() );
newField->SetNumberOfComponents(1);
newField->SetNumberOfTuples(*vlen);

for(unsigned i=0; i<*vlen; i++)
newField->InsertValue(i, vect[i]);

dataSet->GetFieldData()->AddArray(newField);

newField->Delete();
return;
}

/**
Writes float field data values.
@param[in] vect Data array.
@param[in] vlen length of vect.
@param[in] name Variable name to be written to metadata.
@param[in] len Length of variable name.
*/
void vtkwriteffd(float *vect, int *vlen, char *name, int *len){
string tag(name, *len);
vtkFloatArray *newField = vtkFloatArray::New();
newField->SetName( tag.c_str() );
newField->SetNumberOfComponents(1);
newField->SetNumberOfTuples(*vlen);

for(unsigned i=0; i<*vlen; i++)
newField->InsertValue(i, vect[i]);

dataSet->GetFieldData()->AddArray(newField);

newField->Delete();
return;
}

/**
Writes double field data values.
@param[in] vect Data array.
@param[in] vlen length of vect.
@param[in] name Variable name to be written to metadata.
@param[in] len Length of variable name.
*/
void vtkwritedfd(double *vect, int *vlen, char *name, int *len){
string tag(name, *len);
vtkDoubleArray *newField = vtkDoubleArray::New();
newField->SetName( tag.c_str() );
newField->SetNumberOfComponents(1);
newField->SetNumberOfTuples(*vlen);

for(unsigned i=0; i<*vlen; i++)
newField->InsertValue(i, vect[i]);

dataSet->GetFieldData()->AddArray(newField);

newField->Delete();
return;
}


/**
Writes double field data values.
@param[in] vect Data array.
@param[in] vlen length of vect.
@param[in] name Variable name to be written to metadata.
@param[in] len Length of variable name.
*/
void vtkwritecfd(char *vect, int *vlen, char *name, int *len){
string tag(name, *len);
vtkUnsignedCharArray *newField = vtkUnsignedCharArray::New();
newField->SetName( tag.c_str() );
newField->SetNumberOfComponents(1);
newField->SetNumberOfTuples(*vlen);

for(unsigned i=0; i<*vlen; i++)
newField->InsertValue(i, vect[i]);

dataSet->GetFieldData()->AddArray(newField);

newField->Delete();
return;
}


/**
Finish writing and close vtk file (serial).
*/
Expand Down
6 changes: 6 additions & 0 deletions schemas/fluidity_options.rnc
Original file line number Diff line number Diff line change
Expand Up @@ -323,6 +323,12 @@ start =
element memory_diagnostics {
empty
}?
}?,
## include time as VTK FieldData with the name "TimeValue".
## Paraview versions >= 5.5 will automatically use this to
## set the time level.
element time_as_metadata {
empty
}?
},
## Particle options. Particles, like Lagrangian detectors, move with a specified velocity field.
Expand Down
8 changes: 8 additions & 0 deletions schemas/fluidity_options.rng
Original file line number Diff line number Diff line change
Expand Up @@ -399,6 +399,14 @@ NOTE: Requires -v2</a:documentation>
</optional>
</element>
</optional>
<optional>
<element name="time_as_metadata">
<a:documentation>include time as VTK FieldData with the name "TimeValue".
Paraview versions &gt;= 5.5 will automatically use this to
set the time level.</a:documentation>
<empty/>
</element>
</optional>
</element>
<optional>
<element name="particles">
Expand Down