Skip to content
Merged
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
1 change: 1 addition & 0 deletions .github/workflows/ubuntu.yml
Original file line number Diff line number Diff line change
Expand Up @@ -88,6 +88,7 @@ jobs:

- name: Run python tests
env:
SEMBA_FDTD_ENABLE_MPI: ${{ matrix.mpi }}
SEMBA_FDTD_ENABLE_MTLN: ${{ matrix.mtln }}
SEMBA_FDTD_ENABLE_HDF: ${{ matrix.hdf }}
run: python -m pytest test/
Expand Down
1 change: 1 addition & 0 deletions .github/workflows/windows.yml
Original file line number Diff line number Diff line change
Expand Up @@ -69,6 +69,7 @@ jobs:

- name: Run python tests (except codemodel)
env:
SEMBA_FDTD_ENABLE_MPI: ${{ matrix.mpi }}
SEMBA_FDTD_ENABLE_MTLN: ${{ matrix.mtln }}
SEMBA_FDTD_ENABLE_HDF: ${{ matrix.hdf }}
run: python -m pytest -m 'not codemodel' test/
Expand Down
54 changes: 49 additions & 5 deletions src_main_pub/timestepping.F90
Original file line number Diff line number Diff line change
@@ -1,3 +1,44 @@
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! MIT License
!
! Copyright (c) 2023 University of Granada
!
! Permission is hereby granted, free of charge, to any person obtaining a copy
! of this software and associated documentation files (the "Software"), to deal
! in the Software without restriction, including without limitation the rights
! to use, copy, modify, merge, publish, distribute, sublicense, and/or sell
! copies of the Software, and to permit persons to whom the Software is
! furnished to do so, subject to the following conditions:
!
! The above copyright notice and this permission notice shall be included in all
! copies or substantial portions of the Software.
!
! THE SOFTWARE IS PROVIDED "AS IS", WITHOUT WARRANTY OF ANY KIND, EXPRESS OR
! IMPLIED, INCLUDING BUT NOT LIMITED TO THE WARRANTIES OF MERCHANTABILITY,
! FITNESS FOR A PARTICULAR PURPOSE AND NONINFRINGEMENT. IN NO EVENT SHALL THE
! AUTHORS OR COPYRIGHT HOLDERS BE LIABLE FOR ANY CLAIM, DAMAGES OR OTHER
! LIABILITY, WHETHER IN AN ACTION OF CONTRACT, TORT OR OTHERWISE, ARISING FROM,
! OUT OF OR IN CONNECTION WITH THE SOFTWARE OR THE USE OR OTHER DEALINGS IN THE
! SOFTWARE.
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!


!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!
! SEMBA_FDTD sOLVER MODULE
! Creation date Date : April, 8, 2010
!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!!

!__________________________________________________________________________________________________
!******************************** REVISAR PARA PGI (CRAY) *****************************************
!---> AdvanceMultiportE
!---> AdvanceAnisMultiportE
!---> AdvanceMultiportH
!---> AdvanceAnisMultiportH
!---> dfUpdateE
!---> dfUpdateH
!---> MinusCloneMagneticPMC
!________________________________________________________________________________________

module Solver

use fdetypes
Expand Down Expand Up @@ -2782,8 +2823,11 @@ subroutine fillMtag(sgg,sggMiEx, sggMiEy, sggMiEz, sggMiHx, sggMiHy, sggMiHz,sgg
integer(kind = 4) :: i, j, k
integer(kind = INTEGERSIZEOFMEDIAMATRICES) :: medio1,medio2,medio3,medio4,medio5
logical :: mediois1,mediois2,mediois3,mediois4


integer, dimension(3) :: lbx, lby, lbz
lbx = lbound(tag_numbers%face%x)
lby = lbound(tag_numbers%face%y)
lbz = lbound(tag_numbers%face%z)

mediois3=.true.; mediois4=.true.
#ifdef CompileWithOpenMP
!$OMP PARALLEL DO DEFAULT(SHARED) private (i,j,k,medio1,medio2,medio3,medio4,medio5,mediois1,mediois2,mediois3,mediois4)
Expand All @@ -2801,7 +2845,7 @@ subroutine fillMtag(sgg,sggMiEx, sggMiEy, sggMiEz, sggMiHx, sggMiHy, sggMiHz,sgg
mediois3= .true. !.not.((medio5==1).and.(((sggMiHx(i-1,j,k)/=1).or.(sggMiHx(i+1,j,k)/=1)))) !esta condicion en realidad no detecta alabeos de una celda que siendo slots son acoples de un agujerito solo en el peor de los casos
if ((mediois1.or.mediois2).and.(mediois3)) then
!solo lo hace con celdas de vacio porque en particular el mismo medio sgbc con diferentes orientaciones tiene distintos indices de medio y lo activaria erroneamente si lo hago para todos los medios
tag_numbers%face%x(i,j,k)=-ibset(iabs(tag_numbers%face%x(i,j,k)),3)
tag_numbers%face%x(i+lbx(1)-1,j+lbx(2)-1,k+lbx(3)-1)=-ibset(iabs(tag_numbers%face%x(i+lbx(1)-1,j+lbx(2)-1,k+lbx(3)-1)),3)
!ojo no cambiar: interacciona con observation tags 141020 !151020 a efectos de mapvtk el signo importa
endif
End do
Expand All @@ -2823,7 +2867,7 @@ subroutine fillMtag(sgg,sggMiEx, sggMiEy, sggMiEz, sggMiHx, sggMiHy, sggMiHz,sgg
mediois2= (medio5==1).and.(medio3/=1).and.(medio4/=1).and.(medio1==1).and.(medio2==1)
mediois3= .true. !.not.((medio5==1).and.(((sggMiHy(i,j-1,k)/=1).or.(sggMiHy(i,j+1,k)/=1))))
if ((mediois1.or.mediois2).and.(mediois3)) then
tag_numbers%face%y(i,j,k)=-ibset(iabs(tag_numbers%face%y(i,j,k)),4)
tag_numbers%face%y(i+lby(1)-1,j+lby(2)-1,k+lby(3)-1)=-ibset(iabs(tag_numbers%face%y(i+lby(1)-1,j+lby(2)-1,k+lby(3)-1)),4)
endif
End do
End do
Expand All @@ -2844,7 +2888,7 @@ subroutine fillMtag(sgg,sggMiEx, sggMiEy, sggMiEz, sggMiHx, sggMiHy, sggMiHz,sgg
mediois2= (medio5==1).and.(medio3/=1).and.(medio4/=1).and.(medio1==1).and.(medio2==1)
mediois3= .true. !.not.((medio5==1).and.(((sggMiHz(i,j,k-1)/=1).or.(sggMiHz(i,j,k+1)/=1))))
if ((mediois1.or.mediois2).and.(mediois3)) then
tag_numbers%face%z(i,j,k)=-ibset(iabs(tag_numbers%face%z(i,j,k)),5)
tag_numbers%face%z(i+lbz(1)-1,j+lbz(2)-1,k+lbz(3)-1)=-ibset(iabs(tag_numbers%face%z(i+lbz(1)-1,j+lbz(2)-1,k+lbz(3)-1)),5)
endif
End do
End do
Expand Down
36 changes: 25 additions & 11 deletions src_pyWrapper/pyWrapper.py
Original file line number Diff line number Diff line change
Expand Up @@ -24,7 +24,8 @@ class Probe():
+ BULK_CURRENT_PROBE_TAGS \
+ POINT_PROBE_TAGS \
+ FAR_FIELD_TAG \
+ MOVIE_TAGS
+ MOVIE_TAGS

def __init__(self, probe_filename):
if isinstance(probe_filename, os.PathLike):
self.filename = probe_filename.as_posix()
Expand Down Expand Up @@ -92,7 +93,7 @@ def __init__(self, probe_filename):
self.field, self.direction = Probe._getFieldAndDirection(tag)
self.cell_init, self.cell_end = \
Probe._positionStrToTwoCells(position_str)

if self.domainType == 'time':
self.data = self.data.rename(columns={
't': 'time',
Expand Down Expand Up @@ -179,17 +180,32 @@ def _getFieldAndDirection(tag: str):


class FDTD():
def __init__(self, input_filename, path_to_exe=None, flags=[], run_in_folder=None):
def __init__(self, input_filename, path_to_exe=None,
flags=None, run_in_folder=None, mpi_command=None):

self._setFilename(input_filename)

if path_to_exe is None:
self.path_to_exe = os.path.join(
os.getcwd(), DEFAULT_SEMBA_FDTD_PATH)
semba_exe = \
os.path.join(os.getcwd(), DEFAULT_SEMBA_FDTD_PATH)
else:
semba_exe = path_to_exe
assert os.path.isfile(semba_exe)

if mpi_command is None:
mpi_command_parts = []
else:
self.path_to_exe = path_to_exe
assert os.path.isfile(self.path_to_exe)
mpi_command_parts = mpi_command.split()

if flags is None:
flags = []
elif isinstance(flags, str):
flags = flags.split()

case_name = self.getCaseName() + ".json"
self.run_command = \
mpi_command_parts + [semba_exe] + ["-i", case_name] + flags

self.flags = flags
self._hasRun = False

if run_in_folder != None:
Expand Down Expand Up @@ -263,9 +279,7 @@ def run(self):
json.dump(self._input, open(self._filename, 'w'))

os.chdir(self.getFolder())
case_name = self.getCaseName() + ".json"
self.output = subprocess.run(
[self.path_to_exe, "-i", case_name]+self.flags)
self.output = subprocess.run(self.run_command)

self._hasRun = True
assert self.hasFinishedSuccessfully()
Expand Down
148 changes: 93 additions & 55 deletions test/pyWrapper/test_integration.py
Original file line number Diff line number Diff line change
Expand Up @@ -39,6 +39,19 @@ def test_towel_hanger_case_creates_output_probes(tmp_path):
assert countLinesInFile(probe_mid[0]) == 3
assert countLinesInFile(probe_end[0]) == 3

@no_mpi_skip
def test_airplane_case_with_mpi(tmp_path):
fn = CASES_FOLDER + 'airplane/airplane.fdtd.json'
solver = FDTD(fn,
path_to_exe=SEMBA_EXE,
run_in_folder=tmp_path,
flags=['-mapvtk'],
mpi_command='mpirun -np 2')
solver.run()

vtkmapfile = solver.getVTKMap()
assert os.path.isfile(vtkmapfile)


def test_sphere_case_with_far_field_probe_launches(tmp_path):
fn = CASES_FOLDER + 'sphere/sphere.fdtd.json'
Expand Down Expand Up @@ -66,29 +79,34 @@ def test_tagnumbers_3_surfaces(tmp_path):
solver['general']['numberOfSteps'] = 1

solver.run()

vtkmapfile = solver.getVTKMap()
assert os.path.isfile(vtkmapfile)

face_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'tagnumber')
face_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='tagnumber')
assert face_tag_dict[64] == 4
assert face_tag_dict[128] == 4
assert face_tag_dict[192] == 4

line_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'tagnumber')
line_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='tagnumber')
assert line_tag_dict[64] == 8
assert line_tag_dict[128] == 4
assert line_tag_dict[192] == 4

face_media_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'mediatype')
assert face_media_dict[0] == 4 #PEC surface
assert face_media_dict[304] == 4 #SGBC surface
assert face_media_dict[305] == 4 #SGBC surface

line_media_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'mediatype')
assert line_media_dict[0.5] == 8 #PEC line
assert line_media_dict[3.5] == 8 #SGBC line


face_media_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='mediatype')
assert face_media_dict[0] == 4 # PEC surface
assert face_media_dict[304] == 4 # SGBC surface
assert face_media_dict[305] == 4 # SGBC surface

line_media_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='mediatype')
assert line_media_dict[0.5] == 8 # PEC line
assert line_media_dict[3.5] == 8 # SGBC line


def test_tagnumbers_1_volume(tmp_path):
fn = CASES_FOLDER + 'tagNumber_mediaType/pec_volume.fdtd.json'
solver = FDTD(input_filename=fn, path_to_exe=SEMBA_EXE,
Expand All @@ -100,18 +118,23 @@ def test_tagnumbers_1_volume(tmp_path):
vtkmapfile = solver.getVTKMap()
assert os.path.isfile(vtkmapfile)

face_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'tagnumber')
face_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='tagnumber')
assert face_tag_dict[64] == 36

line_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'tagnumber')

line_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='tagnumber')
assert len(line_tag_dict) == 0

face_media_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'mediatype')
assert face_media_dict[0] == 36 #PEC surface

line_media_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'mediatype')

face_media_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='mediatype')
assert face_media_dict[0] == 36 # PEC surface

line_media_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='mediatype')
assert len(line_media_dict) == 0


def test_tagnumbers_2_volumes(tmp_path):
fn = CASES_FOLDER + 'tagNumber_mediaType/pec_volumes.fdtd.json'
solver = FDTD(input_filename=fn, path_to_exe=SEMBA_EXE,
Expand All @@ -123,19 +146,24 @@ def test_tagnumbers_2_volumes(tmp_path):
vtkmapfile = solver.getVTKMap()
assert os.path.isfile(vtkmapfile)

face_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'tagnumber')
face_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='tagnumber')
assert face_tag_dict[64] == 36
assert face_tag_dict[128] == 36

line_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'tagnumber')

line_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='tagnumber')
assert len(line_tag_dict) == 0

face_media_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'mediatype')
assert face_media_dict[0] == 72 #PEC surface

line_media_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'mediatype')

face_media_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='mediatype')
assert face_media_dict[0] == 72 # PEC surface

line_media_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='mediatype')
assert len(line_media_dict) == 0


def test_tagnumbers_1_line(tmp_path):
fn = CASES_FOLDER + 'tagNumber_mediaType/pec_line.fdtd.json'
solver = FDTD(input_filename=fn, path_to_exe=SEMBA_EXE,
Expand All @@ -147,18 +175,23 @@ def test_tagnumbers_1_line(tmp_path):
vtkmapfile = solver.getVTKMap()
assert os.path.isfile(vtkmapfile)

face_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'tagnumber')
face_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='tagnumber')
assert len(face_tag_dict) == 0

line_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'tagnumber')

line_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='tagnumber')
assert line_tag_dict[64] == 2

face_media_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'mediatype')

face_media_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='mediatype')
assert len(face_media_dict) == 0

line_media_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'mediatype')
assert line_media_dict[0.5] == 2 #PEC line


line_media_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='mediatype')
assert line_media_dict[0.5] == 2 # PEC line


def test_tagnumbers_volume_and_surfaces(tmp_path):
fn = CASES_FOLDER + 'tagNumber_mediaType/volume_and_surfaces.fdtd.json'
solver = FDTD(input_filename=fn, path_to_exe=SEMBA_EXE,
Expand All @@ -170,26 +203,31 @@ def test_tagnumbers_volume_and_surfaces(tmp_path):
vtkmapfile = solver.getVTKMap()
assert os.path.isfile(vtkmapfile)

face_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'tagnumber')
face_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='tagnumber')
assert face_tag_dict[64] == 6
assert face_tag_dict[128] == 1
assert face_tag_dict[192] == 1

line_tag_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'tagnumber')

line_tag_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='tagnumber')
assert line_tag_dict[64] == 1
assert line_tag_dict[128] == 4
assert line_tag_dict[192] == 3

face_media_dict = createPropertyDictionary(vtkmapfile, celltype = 9, property = 'mediatype')
assert face_media_dict[-1] == 1 #PEC surface
assert face_media_dict[0] == 6 #PEC surface
assert face_media_dict[305] == 1 #SGBC surface

line_media_dict = createPropertyDictionary(vtkmapfile, celltype = 3, property = 'mediatype')
assert line_media_dict[-0.5] == 4 #PMC line
assert line_media_dict[0.5] == 1 #PEC line
assert line_media_dict[3.5] == 3 #SGBC line


face_media_dict = createPropertyDictionary(
vtkmapfile, celltype=9, property='mediatype')
assert face_media_dict[-1] == 1 # PEC surface
assert face_media_dict[0] == 6 # PEC surface
assert face_media_dict[305] == 1 # SGBC surface

line_media_dict = createPropertyDictionary(
vtkmapfile, celltype=3, property='mediatype')
assert line_media_dict[-0.5] == 4 # PMC line
assert line_media_dict[0.5] == 1 # PEC line
assert line_media_dict[3.5] == 3 # SGBC line


def test_tagnumbers_count_bug(tmp_path):
fn = CASES_FOLDER + 'tagNumber_mediaType/count_bug.fdtd.json'
solver = FDTD(input_filename=fn, path_to_exe=SEMBA_EXE,
Expand All @@ -201,10 +239,10 @@ def test_tagnumbers_count_bug(tmp_path):
solver["materialAssociations"][1]["materialId"] = 1
solver["materialAssociations"][2]["materialId"] = 3
solver.cleanUp()
solver.run()
solver.run()

solver["materialAssociations"][0]["materialId"] = 3
solver["materialAssociations"][1]["materialId"] = 3
solver["materialAssociations"][2]["materialId"] = 1
solver.cleanUp()
solver.run()
solver.run()
Loading